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/**
166 * @brief Check if a string represents a valid number.
167 *
168 * This function checks whether the provided string 'str' represents a valid number.
169 * The function returns true if the entire string can be parsed as a number (integer or floating-point); otherwise, it
170 * returns false.
171 *
172 * @param str The string to be checked.
173 * @return bool Returns true if the string 'str' represents a valid number; otherwise, returns false.
174 */
175bool isNumber(const std::string &str)
176{
177 bool seen_digit = false;
178 bool seen_dot = false;
179 bool seen_e = false;
180 bool after_e = false;
181 bool sign_allowed = true;
182
183 for (size_t i = 0; i < str.size(); ++i) {
184 char c = str[i];
185
186 if (std::isdigit(c)) {
187 seen_digit = true;
188 sign_allowed = false;
189 } else if ((c == '+' || c == '-') && sign_allowed) {
190 // Sign allowed at the beginning or right after 'e'/'E'
191 sign_allowed = false;
192 } else if (c == '.' && !seen_dot && !after_e) {
193 seen_dot = true;
194 sign_allowed = false;
195 } else if ((c == 'e' || c == 'E') && seen_digit && !seen_e) {
196 seen_e = true;
197 after_e = true;
198 sign_allowed = true; // allow sign immediately after 'e'
199 seen_digit = false; // reset: we now expect digits after e
200 } else {
201 return false;
202 }
203 }
204
205 return seen_digit;
206}
207
208/**
209 * @brief Configure a RooRealVar based on information from a JSONNode.
210 *
211 * This function configures the provided RooRealVar 'v' based on the information provided in the JSONNode 'p'.
212 * The JSONNode 'p' contains information about various properties of the RooRealVar, such as its value, error, number of
213 * bins, etc. The function reads these properties from the JSONNode and sets the corresponding properties of the
214 * RooRealVar accordingly.
215 *
216 * @param domains The reference to the RooFit::JSONIO::Detail::Domains containing domain information for variables (not
217 * used in this function).
218 * @param p The JSONNode containing information about the properties of the RooRealVar 'v'.
219 * @param v The reference to the RooRealVar to be configured.
220 * @return void
221 */
223{
224 if (!p.has_child("name")) {
225 RooJSONFactoryWSTool::error("cannot instantiate variable without \"name\"!");
226 }
227 if (auto n = p.find("value"))
228 v.setVal(n->val_double());
229 domains.writeVariable(v);
230 if (auto n = p.find("nbins"))
231 v.setBins(n->val_int());
232 if (auto n = p.find("relErr"))
233 v.setError(v.getVal() * n->val_double());
234 if (auto n = p.find("err"))
235 v.setError(n->val_double());
236 if (auto n = p.find("const")) {
237 v.setConstant(n->val_bool());
238 } else {
239 v.setConstant(false);
240 }
241}
242
244{
245 auto paramPointsNode = rootNode.find("parameter_points");
246 if (!paramPointsNode)
247 return nullptr;
248 auto out = RooJSONFactoryWSTool::findNamedChild(*paramPointsNode, "default_values");
249 if (out == nullptr)
250 return nullptr;
251 return &((*out)["parameters"]);
252}
253
254std::string genPrefix(const JSONNode &p, bool trailing_underscore)
255{
256 std::string prefix;
257 if (!p.is_map())
258 return prefix;
259 if (auto node = p.find("namespaces")) {
260 for (const auto &ns : node->children()) {
261 if (!prefix.empty())
262 prefix += "_";
263 prefix += ns.val();
264 }
265 }
266 if (trailing_underscore && !prefix.empty())
267 prefix += "_";
268 return prefix;
269}
270
271// helpers for serializing / deserializing binned datasets
272void genIndicesHelper(std::vector<std::vector<int>> &combinations, std::vector<int> &curr_comb,
273 const std::vector<int> &vars_numbins, size_t curridx)
274{
275 if (curridx == vars_numbins.size()) {
276 // we have filled a combination. Copy it.
277 combinations.emplace_back(curr_comb);
278 } else {
279 for (int i = 0; i < vars_numbins[curridx]; ++i) {
280 curr_comb[curridx] = i;
282 }
283 }
284}
285
286/**
287 * @brief Import attributes from a JSONNode into a RooAbsArg.
288 *
289 * This function imports attributes, represented by the provided JSONNode 'node', into the provided RooAbsArg 'arg'.
290 * The attributes are read from the JSONNode and applied to the RooAbsArg.
291 *
292 * @param arg The pointer to the RooAbsArg to which the attributes will be imported.
293 * @param node The JSONNode containing information about the attributes to be imported.
294 * @return void
295 */
296void importAttributes(RooAbsArg *arg, JSONNode const &node)
297{
298 if (auto seq = node.find("dict")) {
299 for (const auto &attr : seq->children()) {
300 arg->setStringAttribute(attr.key().c_str(), attr.val().c_str());
301 }
302 }
303 if (auto seq = node.find("tags")) {
304 for (const auto &attr : seq->children()) {
305 arg->setAttribute(attr.val().c_str());
306 }
307 }
308}
309
310void addIfPresent(RooArgSet &out, RooArgSet const *args)
311{
312 if (args) {
313 out.add(*args, true);
314 }
315}
316
319{
320 for (TObject *obj : workspace.allGenericObjects()) {
321 auto const *mc = dynamic_cast<RooFit::ModelConfig const *>(obj);
322 if (!mc) {
323 continue;
324 }
325
326 addIfPresent(candidates, mc->GetParametersOfInterest());
327 addIfPresent(candidates, mc->GetNuisanceParameters());
328
329 addIfPresent(excluded, mc->GetObservables());
330 addIfPresent(excluded, mc->GetGlobalObservables());
331 addIfPresent(excluded, mc->GetConditionalObservables());
332 }
333}
334
335void collectParameterStepWidthCandidatesFromPdfs(std::vector<RooAbsPdf *> const &pdfs,
336 std::vector<RooAbsData *> const &data, RooArgSet &candidates,
338{
339 for (RooAbsPdf const *pdf : pdfs) {
340 RooArgSet observables;
341 for (RooAbsData const *dataset : data) {
342 std::unique_ptr<RooArgSet> pdfObs{pdf->getObservables(*dataset->get())};
343 observables.add(*pdfObs, true);
344 }
345
346 if (observables.empty()) {
347 continue;
348 }
349
350 RooArgSet params;
351 pdf->getParameters(&observables, params);
352 candidates.add(params, true);
353 excluded.add(observables, true);
354 }
355}
356
357void exportParameterStepWidths(RooWorkspace const &workspace, std::vector<RooAbsPdf *> const &pdfs,
358 std::vector<RooAbsData *> const &data, JSONNode &rootnode)
359{
362
365
366 candidates.sort();
367
369 for (RooAbsArg *arg : candidates) {
370 if (excluded.find(*arg)) {
371 continue;
372 }
373
374 auto *var = dynamic_cast<RooRealVar *>(arg);
375 if (!var || !var->hasError()) {
376 continue;
377 }
378
380 parameterStepWidthsNode = &rootnode["misc"]["minimization"]["parameter_stepwidths"].set_seq();
381 }
382
384 stepWidthNode["step_width"] << var->getError();
385 }
386}
387
388void importParameterStepWidths(RooWorkspace &workspace, JSONNode const &rootnode)
389{
390 auto const *parameterStepWidthsNode = rootnode.find("misc", "minimization", "parameter_stepwidths");
392 return;
393 }
394 if (!parameterStepWidthsNode->is_seq()) {
395 RooJSONFactoryWSTool::warning("RooFitHS3: misc.minimization.parameter_stepwidths is not a sequence, skipping.");
396 return;
397 }
398
399 for (JSONNode const &stepWidthNode : parameterStepWidthsNode->children()) {
400 if (!stepWidthNode.is_map() || !stepWidthNode.has_child("name") || !stepWidthNode.has_child("step_width")) {
401 RooJSONFactoryWSTool::warning("RooFitHS3: skipping malformed parameter_stepwidths entry.");
402 continue;
403 }
404
405 const std::string name = RooJSONFactoryWSTool::name(stepWidthNode);
406 RooAbsArg *arg = workspace.arg(name);
407 auto *var = dynamic_cast<RooRealVar *>(arg);
408 if (!var) {
410 "RooFitHS3: skipping parameter_stepwidths entry for unknown or non-real variable '" + name + "'.");
411 continue;
412 }
413
414 var->setError(stepWidthNode.find("step_width")->val_double());
415 }
416}
417
418// RooWSFactoryTool expression handling
419std::string generate(const RooFit::JSONIO::ImportExpression &ex, const JSONNode &p, RooJSONFactoryWSTool *tool)
420{
421 std::stringstream expression;
422 std::string classname(ex.tclass->GetName());
423 size_t colon = classname.find_last_of(':');
424 expression << (colon < classname.size() ? classname.substr(colon + 1) : classname);
425 bool first = true;
426 const auto &name = RooJSONFactoryWSTool::name(p);
427 for (auto k : ex.arguments) {
428 expression << (first ? "::" + name + "(" : ",");
429 first = false;
430 if (k == "true" || k == "false") {
431 expression << (k == "true" ? "1" : "0");
432 } else if (!p.has_child(k)) {
433 std::stringstream errMsg;
434 errMsg << "node '" << name << "' is missing key '" << k << "'";
436 } else if (p[k].is_seq()) {
437 bool firstInner = true;
438 expression << "{";
439 for (RooAbsArg *arg : tool->requestArgList<RooAbsReal>(p, k)) {
440 expression << (firstInner ? "" : ",") << arg->GetName();
441 firstInner = false;
442 }
443 expression << "}";
444 } else {
445 tool->requestArg<RooAbsReal>(p, p[k].key());
446 expression << p[k].val();
447 }
448 }
449 expression << ")";
450 return expression.str();
451}
452
453/**
454 * @brief Generate bin indices for a set of RooRealVars.
455 *
456 * This function generates all possible combinations of bin indices for the provided RooArgSet 'vars' containing
457 * RooRealVars. Each bin index represents a possible bin selection for the corresponding RooRealVar. The bin indices are
458 * stored in a vector of vectors, where each inner vector represents a combination of bin indices for all RooRealVars.
459 *
460 * @param vars The RooArgSet containing the RooRealVars for which bin indices will be generated.
461 * @return std::vector<std::vector<int>> A vector of vectors containing all possible combinations of bin indices.
462 */
463std::vector<std::vector<int>> generateBinIndices(const RooArgSet &vars)
464{
465 std::vector<std::vector<int>> combinations;
466 std::vector<int> vars_numbins;
467 vars_numbins.reserve(vars.size());
468 for (const auto *absv : static_range_cast<RooRealVar *>(vars)) {
469 vars_numbins.push_back(absv->getBins());
470 }
471 std::vector<int> curr_comb(vars.size());
473 return combinations;
474}
475
476template <typename... Keys_t>
477JSONNode const *findRooFitInternal(JSONNode const &node, Keys_t const &...keys)
478{
479 return node.find("misc", "ROOT_internal", keys...);
480}
481
482/**
483 * @brief Check if a RooAbsArg is a literal constant variable.
484 *
485 * This function checks whether the provided RooAbsArg 'arg' is a literal constant variable.
486 * A literal constant variable is a RooConstVar with a numeric value as a name.
487 *
488 * @param arg The reference to the RooAbsArg to be checked.
489 * @return bool Returns true if 'arg' is a literal constant variable; otherwise, returns false.
490 */
491bool isLiteralConstVar(RooAbsArg const &arg)
492{
493 bool isRooConstVar = dynamic_cast<RooConstVar const *>(&arg);
494 return isRooConstVar && isNumber(arg.GetName());
495}
496
497/**
498 * @brief Export attributes of a RooAbsArg to a JSONNode.
499 *
500 * This function exports the attributes of the provided RooAbsArg 'arg' to the JSONNode 'rootnode'.
501 *
502 * @param arg The pointer to the RooAbsArg from which attributes will be exported.
503 * @param rootnode The JSONNode to which the attributes will be exported.
504 * @return void
505 */
506void exportAttributes(const RooAbsArg *arg, JSONNode &rootnode)
507{
508 // If this RooConst is a literal number, we don't need to export the attributes.
509 if (isLiteralConstVar(*arg)) {
510 return;
511 }
512
513 JSONNode *node = nullptr;
514
515 auto initializeNode = [&]() {
516 if (node)
517 return;
518
519 node = &RooJSONFactoryWSTool::getRooFitInternal(rootnode, "attributes").set_map()[arg->GetName()].set_map();
520 };
521
522 // RooConstVars are not a thing in HS3, and also for RooFit they are not
523 // that important: they are just constants. So we don't need to remember
524 // any information about them.
525 if (dynamic_cast<RooConstVar const *>(arg)) {
526 return;
527 }
528
529 // export all string attributes of an object
530 if (!arg->stringAttributes().empty()) {
531 for (const auto &it : arg->stringAttributes()) {
532 // Skip some RooFit internals
533 if (it.first == "factory_tag" || it.first == "PROD_TERM_TYPE")
534 continue;
536 (*node)["dict"].set_map()[it.first] << it.second;
537 }
538 }
539 if (!arg->attributes().empty()) {
540 for (auto const &attr : arg->attributes()) {
541 // Skip some RooFit internals
542 if (attr == "SnapShot_ExtRefClone" || attr == "RooRealConstant_Factory_Object")
543 continue;
545 (*node)["tags"].set_seq().append_child() << attr;
546 }
547 }
548}
549
550/**
551 * @brief Create several observables in the workspace.
552 *
553 * This function obtains a list of observables from the provided
554 * RooWorkspace 'ws' based on their names given in the 'axes" field of
555 * the JSONNode 'node'. The observables are added to the RooArgSet
556 * 'out'.
557 *
558 * @param ws The RooWorkspace in which the observables will be created.
559 * @param node The JSONNode containing information about the observables to be created.
560 * @param out The RooAbsCollection to which the created observables will be added.
561 * @return void
562 */
563void getObservables(RooWorkspace const &ws, const JSONNode &node, RooAbsCollection &out)
564{
565 for (const auto &p : node["axes"].children()) {
566 std::string name(RooJSONFactoryWSTool::name(p));
567 if (ws.var(name)) {
568 out.add(*ws.var(name));
569 } else {
570 std::stringstream errMsg;
571 errMsg << "The observable \"" << name << "\" could not be found in the workspace!";
573 }
574 }
575}
576
577/**
578 * @brief Import data from the JSONNode into the workspace.
579 *
580 * This function imports data, represented by the provided JSONNode 'p', into the workspace represented by the provided
581 * RooWorkspace. The data information is read from the JSONNode and added to the workspace.
582 *
583 * @param p The JSONNode representing the data to be imported.
584 * @param workspace The RooWorkspace to which the data will be imported.
585 * @return std::unique_ptr<RooAbsData> A unique pointer to the RooAbsData object representing the imported data.
586 * The caller is responsible for managing the memory of the returned object.
587 */
588std::unique_ptr<RooAbsData> loadData(const JSONNode &p, RooWorkspace &workspace)
589{
590 std::string name(RooJSONFactoryWSTool::name(p));
591
593
594 std::string const &type = p["type"].val();
595 if (type == "binned") {
596 // binned
598 } else if (type == "unbinned") {
599 // unbinned
600 RooArgList varlist;
601 getObservables(workspace, p, varlist);
602 RooArgSet vars(varlist);
603 auto data = std::make_unique<RooDataSet>(name, name, vars, RooFit::WeightVar());
604 auto &coords = p["entries"];
605 if (!coords.is_seq()) {
606 RooJSONFactoryWSTool::error("key 'entries' is not a list!");
607 }
608 std::vector<double> weightVals;
609 if (p.has_child("weights")) {
610 auto &weights = p["weights"];
611 if (coords.num_children() != weights.num_children()) {
612 RooJSONFactoryWSTool::error("inconsistent number of entries and weights!");
613 }
614 for (auto const &weight : weights.children()) {
615 weightVals.push_back(weight.val_double());
616 }
617 }
618 std::size_t i = 0;
619 for (auto const &point : coords.children()) {
620 if (!point.is_seq()) {
621 std::stringstream errMsg;
622 errMsg << "coordinate point '" << i << "' is not a list!";
624 }
625 if (point.num_children() != varlist.size()) {
626 RooJSONFactoryWSTool::error("inconsistent number of entries and observables!");
627 }
628 std::size_t j = 0;
629 for (auto const &pointj : point.children()) {
630 auto *v = static_cast<RooRealVar *>(varlist.at(j));
631 v->setVal(pointj.val_double());
632 ++j;
633 }
634 if (weightVals.size() > 0) {
635 data->add(vars, weightVals[i]);
636 } else {
637 data->add(vars, 1.);
638 }
639 ++i;
640 }
641 return data;
642 }
643
644 std::stringstream ss;
645 ss << "RooJSONFactoryWSTool() failed to create dataset " << name << std::endl;
647 return nullptr;
648}
649
650/**
651 * @brief Import an analysis from the JSONNode into the workspace.
652 *
653 * This function imports an analysis, represented by the provided JSONNodes 'analysisNode' and 'likelihoodsNode',
654 * into the workspace represented by the provided RooWorkspace. The analysis information is read from the JSONNodes
655 * and added to the workspace as one or more RooFit::ModelConfig objects.
656 *
657 * @param rootnode The root JSONNode representing the entire JSON file.
658 * @param analysisNode The JSONNode representing the analysis to be imported.
659 * @param likelihoodsNode The JSONNode containing information about likelihoods associated with the analysis.
660 * @param domainsNode The JSONNode containing information about domains associated with the analysis.
661 * @param workspace The RooWorkspace to which the analysis will be imported.
662 * @param datasets A vector of unique pointers to RooAbsData objects representing the data associated with the analysis.
663 * @return void
664 */
665void importAnalysis(const JSONNode &rootnode, const JSONNode &analysisNode, const JSONNode &likelihoodsNode,
666 const JSONNode &domainsNode, RooWorkspace &workspace,
667 const std::vector<std::unique_ptr<RooAbsData>> &datasets)
668{
669 // if this is a toplevel pdf, also create a modelConfig for it
671 JSONNode const *mcAuxNode = findRooFitInternal(rootnode, "ModelConfigs", analysisName);
672
673 JSONNode const *mcNameNode = mcAuxNode ? mcAuxNode->find("mcName") : nullptr;
674 std::string mcname = mcNameNode ? mcNameNode->val() : analysisName;
675 if (workspace.obj(mcname))
676 return;
677
678 workspace.import(RooFit::ModelConfig{mcname.c_str(), mcname.c_str()});
679 auto *mc = static_cast<RooFit::ModelConfig *>(workspace.obj(mcname));
680 mc->SetWS(workspace);
681
683 if (!nllNode) {
684 throw std::runtime_error("likelihood node not found!");
685 }
686 if (!nllNode->has_child("distributions")) {
687 throw std::runtime_error("likelihood node has no distributions attached!");
688 }
689 if (!nllNode->has_child("data")) {
690 throw std::runtime_error("likelihood node has no data attached!");
691 }
692 std::vector<std::string> nllDistNames = valsToStringVec((*nllNode)["distributions"]);
694 for (auto &nameNode : (*nllNode)["aux_distributions"].children()) {
695 if (RooAbsArg *extConstraint = workspace.arg(nameNode.val())) {
697 }
698 }
699 RooArgSet observables;
700 for (auto &nameNode : (*nllNode)["data"].children()) {
701 bool found = false;
702 for (const auto &d : datasets) {
703 if (d->GetName() == nameNode.val()) {
704 found = true;
705 observables.add(*d->get(), true);
706 }
707 }
708 if (nameNode.val() != "0" && !found)
709 throw std::runtime_error("dataset '" + nameNode.val() + "' cannot be found!");
710 }
711
712 JSONNode const *pdfNameNode = mcAuxNode ? mcAuxNode->find("pdfName") : nullptr;
713 std::string const pdfName = pdfNameNode ? pdfNameNode->val() : "simPdf";
714
715 RooAbsPdf *pdf = static_cast<RooSimultaneous *>(workspace.pdf(pdfName));
716
717 if (!pdf) {
718 // if there is no simultaneous pdf, we can check whether there is only one pdf in the list
719 if (nllDistNames.size() == 1) {
720 // if so, we can use that one to populate the ModelConfig
721 pdf = workspace.pdf(nllDistNames[0]);
722 } else {
723 // otherwise, we have no choice but to build a simPdf by hand
724 std::string simPdfName = analysisName + "_simPdf";
725 std::string indexCatName = analysisName + "_categoryIndex";
726 RooCategory indexCat{indexCatName.c_str(), indexCatName.c_str()};
727 std::map<std::string, RooAbsPdf *> pdfMap;
728 for (std::size_t i = 0; i < nllDistNames.size(); ++i) {
729 indexCat.defineType(nllDistNames[i], i);
730 pdfMap[nllDistNames[i]] = workspace.pdf(nllDistNames[i]);
731 }
732 RooSimultaneous simPdf{simPdfName.c_str(), simPdfName.c_str(), pdfMap, indexCat};
734 pdf = static_cast<RooSimultaneous *>(workspace.pdf(simPdfName));
735 }
736 }
737
738 mc->SetPdf(*pdf);
739
740 if (!extConstraints.empty())
741 mc->SetExternalConstraints(extConstraints);
742
743 auto readArgSet = [&](std::string const &name) {
744 RooArgSet out;
745 for (auto const &child : analysisNode[name].children()) {
746 out.add(*workspace.arg(child.val()));
747 }
748 return out;
749 };
750
751 mc->SetParametersOfInterest(readArgSet("parameters_of_interest"));
752 mc->SetObservables(observables);
753 RooArgSet pars;
754 pdf->getParameters(&observables, pars);
755
756 // Figure out the set parameters that appear in the main measurement:
757 // getAllConstraints() has the side effect to remove all parameters from
758 // "mainPars" that are not part of any pdf over observables.
759 RooArgSet mainPars{pars};
760 pdf->getAllConstraints(observables, mainPars, /*stripDisconnected*/ true);
761
763 for (auto &domain : analysisNode["domains"].children()) {
765 if (!thisDomain || !thisDomain->has_child("axes"))
766 continue;
767 for (auto &var : (*thisDomain)["axes"].children()) {
768 auto *wsvar = workspace.var(RooJSONFactoryWSTool::name(var));
769 if (wsvar)
770 domainPars.add(*wsvar);
771 }
772 }
773
775 RooArgSet globs;
776 for (const auto &p : pars) {
777 if (mc->GetParametersOfInterest()->find(*p))
778 continue;
779 if (p->isConstant() && !mainPars.find(*p) && domainPars.find(*p)) {
780 globs.add(*p);
781 } else if (domainPars.find(*p)) {
782 nps.add(*p);
783 }
784 }
785
786 mc->SetGlobalObservables(globs);
787 mc->SetNuisanceParameters(nps);
788
789 if (mcAuxNode) {
790 if (auto found = mcAuxNode->find("combined_data_name")) {
791 pdf->setStringAttribute("combined_data_name", found->val().c_str());
792 }
793 }
794
795 if (analysisNode.has_child("init") && workspace.getSnapshot(analysisNode["init"].val().c_str())) {
796 mc->SetSnapshot(*workspace.getSnapshot(analysisNode["init"].val().c_str()));
797 }
798}
799
800void combinePdfs(const JSONNode &rootnode, RooWorkspace &ws)
801{
802 auto *combinedPdfInfoNode = findRooFitInternal(rootnode, "combined_distributions");
803
804 // If there is no info on combining pdfs
805 if (combinedPdfInfoNode == nullptr) {
806 return;
807 }
808
809 for (auto &info : combinedPdfInfoNode->children()) {
810
811 // parse the information
812 std::string combinedName = info.key();
813 std::string indexCatName = info["index_cat"].val();
814 std::vector<std::string> labels = valsToStringVec(info["labels"]);
815 std::vector<int> indices;
816 std::vector<std::string> pdfNames = valsToStringVec(info["distributions"]);
817 for (auto &n : info["indices"].children()) {
818 indices.push_back(n.val_int());
819 }
820
821 RooCategory indexCat{indexCatName.c_str(), indexCatName.c_str()};
822 std::map<std::string, RooAbsPdf *> pdfMap;
823
824 for (std::size_t iChannel = 0; iChannel < labels.size(); ++iChannel) {
825 indexCat.defineType(labels[iChannel], indices[iChannel]);
826 pdfMap[labels[iChannel]] = ws.pdf(pdfNames[iChannel]);
827 }
828
829 RooSimultaneous simPdf{combinedName.c_str(), combinedName.c_str(), pdfMap, indexCat};
831 }
832}
833
834void combineDatasets(const JSONNode &rootnode, std::vector<std::unique_ptr<RooAbsData>> &datasets)
835{
836 auto *combinedDataInfoNode = findRooFitInternal(rootnode, "combined_datasets");
837
838 // If there is no info on combining datasets
839 if (combinedDataInfoNode == nullptr) {
840 return;
841 }
842
843 for (auto &info : combinedDataInfoNode->children()) {
844
845 // parse the information
846 std::string combinedName = info.key();
847 std::string indexCatName = info["index_cat"].val();
848 std::vector<std::string> labels = valsToStringVec(info["labels"]);
849 std::vector<int> indices;
850 for (auto &n : info["indices"].children()) {
851 indices.push_back(n.val_int());
852 }
853 if (indices.size() != labels.size()) {
854 RooJSONFactoryWSTool::error("mismatch in number of indices and labels!");
855 }
856
857 // Create the combined dataset for RooFit
858 std::map<std::string, std::unique_ptr<RooAbsData>> dsMap;
859 RooCategory indexCat{indexCatName.c_str(), indexCatName.c_str()};
860 RooArgSet allVars{indexCat};
861 for (std::size_t iChannel = 0; iChannel < labels.size(); ++iChannel) {
862 auto componentName = combinedName + "_" + labels[iChannel];
863 // We move the found channel data out of the "datasets" vector, such that
864 // the data components don't get imported anymore.
865 std::unique_ptr<RooAbsData> &component = *std::find_if(
866 datasets.begin(), datasets.end(), [&](auto &d) { return d && d->GetName() == componentName; });
867 if (!component)
868 RooJSONFactoryWSTool::error("unable to obtain component matching component name '" + componentName + "'");
869 allVars.add(*component->get(), true);
870 dsMap.insert({labels[iChannel], std::move(component)});
871 indexCat.defineType(labels[iChannel], indices[iChannel]);
872 }
873
874 auto combined = std::make_unique<RooDataSet>(combinedName, combinedName, allVars, RooFit::Import(dsMap),
875 RooFit::Index(indexCat));
876 datasets.emplace_back(std::move(combined));
877 }
878}
879
880template <class T>
881void sortByName(T &coll)
882{
883 std::sort(coll.begin(), coll.end(), [](auto &l, auto &r) { return strcmp(l->GetName(), r->GetName()) < 0; });
884}
885
886} // namespace
887
889
891
893{
894 const size_t old_children = node.num_children();
895 node.set_seq();
896 size_t n = 0;
897 for (RooAbsArg const *arg : coll) {
898 if (n >= nMax)
899 break;
900 if (isLiteralConstVar(*arg)) {
901 node.append_child() << static_cast<RooConstVar const *>(arg)->getVal();
902 } else {
903 node.append_child() << arg->GetName();
904 }
905 ++n;
906 }
907 if (node.num_children() != old_children + coll.size()) {
908 error("unable to stream collection " + std::string(coll.GetName()) + " to " + node.key());
909 }
910}
911
913{
914 const size_t old_children = node.num_children();
915 node.set_seq();
916 size_t n = 0;
917 for (RooAbsArg const *arg : coll) {
918 if (n >= nMax)
919 break;
920 if (isLiteralConstVar(*arg)) {
921 node.append_child() << static_cast<RooConstVar const *>(arg)->getVal();
922 } else {
923 node.append_child() << sanitizeName(arg->GetName());
924 }
925 ++n;
926 }
927 if (node.num_children() != old_children + coll.size()) {
928 error("unable to stream collection " + std::string(coll.GetName()) + " to " + node.key());
929 }
930}
931
933{
935 return node.set_map()[name].set_map();
936 }
937 JSONNode &child = node.set_seq().append_child().set_map();
938 child["name"] << name;
939 return child;
940}
941
942JSONNode const *RooJSONFactoryWSTool::findNamedChild(JSONNode const &node, std::string const &name)
943{
945 if (!node.is_map())
946 return nullptr;
947 return node.find(name);
948 }
949 if (!node.is_seq())
950 return nullptr;
951 for (JSONNode const &child : node.children()) {
952 if (child["name"].val() == name)
953 return &child;
954 }
955
956 return nullptr;
957}
958
959/**
960 * @brief Check if a string is a valid name.
961 *
962 * A valid name should start with a letter or an underscore, followed by letters, digits, or underscores.
963 * Only characters from the ASCII character set are allowed.
964 *
965 * @param str The string to be checked.
966 * @return bool Returns true if the string is a valid name; otherwise, returns false.
967 */
968bool RooJSONFactoryWSTool::isValidName(const std::string &str)
969{
970 // Check if the string is empty or starts with a non-letter/non-underscore character
971 if (str.empty() || !(std::isalpha(str[0]) || str[0] == '_')) {
972 return false;
973 }
974
975 // Check the remaining characters in the string
976 for (char c : str) {
977 // Allow letters, digits, and underscore
978 if (!(std::isalnum(c) || c == '_')) {
979 return false;
980 }
981 }
982
983 // If all characters are valid, the string is a valid name
984 return true;
985}
986
990{
992 std::stringstream ss;
993 ss << "RooJSONFactoryWSTool() name '" << name << "' is not valid!" << std::endl
994 << "Sanitize names by setting RooJSONFactoryWSTool::allowSanitizeNames = True." << std::endl;
997 return false;
998 } else {
1000 }
1001 }
1002 return true;
1003}
1004
1006{
1007 return useListsInsteadOfDicts ? n["name"].val() : n.key();
1008}
1009
1011{
1012 return appendNamedChild(rootNode["parameter_points"], "default_values")["parameters"];
1013}
1014
1015template <>
1016RooRealVar *RooJSONFactoryWSTool::requestImpl<RooRealVar>(const std::string &objname)
1017{
1019 return retval;
1020 if (const auto *vars = getVariablesNode(*_rootnodeInput)) {
1021 if (const auto &node = vars->find(objname)) {
1022 this->importVariable(*node);
1024 return retval;
1025 }
1026 }
1027 return nullptr;
1028}
1029
1030template <>
1031RooAbsPdf *RooJSONFactoryWSTool::requestImpl<RooAbsPdf>(const std::string &objname)
1032{
1034 return retval;
1035 auto it = _distributionsByName.find(objname);
1036 if (it != _distributionsByName.end()) {
1037 this->importFunction(*it->second, true);
1039 return retval;
1040 }
1041 return nullptr;
1042}
1043
1044template <>
1045RooAbsReal *RooJSONFactoryWSTool::requestImpl<RooAbsReal>(const std::string &objname)
1046{
1048 return retval;
1049 if (isNumber(objname))
1052 return pdf;
1054 return var;
1055 auto it = _functionsByName.find(objname);
1056 if (it != _functionsByName.end()) {
1057 this->importFunction(*it->second, true);
1059 return retval;
1060 }
1061 return nullptr;
1062}
1063
1064/**
1065 * @brief Export a variable from the workspace to a JSONNode.
1066 *
1067 * This function exports a variable, represented by the provided RooAbsArg pointer 'v', from the workspace to a
1068 * JSONNode. The variable's information is added to the JSONNode as key-value pairs.
1069 *
1070 * @param v The pointer to the RooAbsArg representing the variable to be exported.
1071 * @param node The JSONNode to which the variable will be exported.
1072 * @return void
1073 */
1075{
1076 auto *cv = dynamic_cast<const RooConstVar *>(v);
1077 auto *rrv = dynamic_cast<const RooRealVar *>(v);
1078 if (!cv && !rrv)
1079 return;
1080
1081 // for RooConstVar, if name and value are the same, we don't need to do anything
1082 if (cv && strcmp(cv->GetName(), TString::Format("%g", cv->getVal()).Data()) == 0) {
1083 return;
1084 }
1085
1086 JSONNode &var = appendNamedChild(node, v->GetName());
1087
1088 if (cv) {
1089 var["value"] << cv->getVal();
1090 var["const"] << true;
1091 } else if (rrv) {
1092 var["value"] << rrv->getVal();
1093 if (rrv->isConstant() && storeConstant) {
1094 var["const"] << rrv->isConstant();
1095 } else {
1096 var["min"] << rrv->getMin();
1097 var["max"] << rrv->getMax();
1098 }
1099 if (rrv->getBins() != 100 && storeBins) {
1100 var["nbins"] << rrv->getBins();
1101 }
1102 _domains->readVariable(*rrv);
1103 }
1104}
1105
1106/**
1107 * @brief Export variables from the workspace to a JSONNode.
1108 *
1109 * This function exports variables, represented by the provided RooArgSet, from the workspace to a JSONNode.
1110 * The variables' information is added to the JSONNode as key-value pairs.
1111 *
1112 * @param allElems The RooArgSet representing the variables to be exported.
1113 * @param n The JSONNode to which the variables will be exported.
1114 * @return void
1115 */
1117{
1118 // export a list of RooRealVar objects
1119 n.set_seq();
1120 for (RooAbsArg *arg : allElems) {
1122 }
1123}
1124
1126 const std::string &formula)
1127{
1128 std::string newname = std::string(original->GetName()) + suffix;
1130 trafo_node["type"] << "generic_function";
1131 trafo_node["expression"] << TString::Format(formula.c_str(), original->GetName()).Data();
1132 this->setAttribute(newname, "roofit_skip"); // this function should not be imported back in
1133 return newname;
1134}
1135
1136/**
1137 * @brief Export an object from the workspace to a JSONNode.
1138 *
1139 * This function exports an object, represented by the provided RooAbsArg, from the workspace to a JSONNode.
1140 * The object's information is added to the JSONNode as key-value pairs.
1141 *
1142 * @param func The RooAbsArg representing the object to be exported.
1143 * @param exportedObjectNames A set of strings containing names of previously exported objects to avoid duplicates.
1144 * This set is updated with the name of the newly exported object.
1145 * @return void
1146 */
1148{
1149 // const std::string name = sanitizeName(func.GetName());
1150 std::string name = func.GetName();
1151
1152 // if this element was already exported, skip
1154 return;
1155
1156 exportedObjectNames.insert(name);
1157
1158 if (auto simPdf = dynamic_cast<RooSimultaneous const *>(&func)) {
1159 // RooSimultaneous is not used in the HS3 standard, we only export the
1160 // dependents and some ROOT internal information.
1162
1163 std::vector<std::string> channelNames;
1164 for (auto const &item : simPdf->indexCat()) {
1165 channelNames.push_back(item.first);
1166 }
1167
1168 auto &infoNode = getRooFitInternal(*_rootnodeOutput, "combined_distributions").set_map();
1169 auto &child = infoNode[simPdf->GetName()].set_map();
1170 child["index_cat"] << simPdf->indexCat().GetName();
1171 exportCategory(simPdf->indexCat(), child);
1172 child["distributions"].set_seq();
1173 for (auto const &item : simPdf->indexCat()) {
1174 child["distributions"].append_child() << simPdf->getPdf(item.first.c_str())->GetName();
1175 }
1176
1177 return;
1178 } else if (dynamic_cast<RooAbsCategory const *>(&func)) {
1179 // categories are created by the respective RooSimultaneous, so we're skipping the export here
1180 return;
1181 } else if (dynamic_cast<RooRealVar const *>(&func) || dynamic_cast<RooConstVar const *>(&func)) {
1182 exportVariable(&func, *_varsNode, true, false);
1183 return;
1184 }
1185
1186 auto &collectionNode = (*_rootnodeOutput)[dynamic_cast<RooAbsPdf const *>(&func) ? "distributions" : "functions"];
1187
1188 auto const &exporters = RooFit::JSONIO::exporters();
1189 auto const &exportKeys = RooFit::JSONIO::exportKeys();
1190
1191 TClass *cl = func.IsA();
1192
1194
1195 auto it = exporters.find(cl);
1196 if (it != exporters.end()) { // check if we have a specific exporter available
1197 for (auto &exp : it->second) {
1198 _serversToExport.clear();
1199 _serversToDelete.clear();
1200 if (!exp->exportObject(this, &func, elem)) {
1201 // The exporter might have messed with the content of the node
1202 // before failing. That's why we clear it and only reset the name.
1203 elem.clear();
1204 elem.set_map();
1206 elem["name"] << name;
1207 }
1208 continue;
1209 }
1210 if (exp->autoExportDependants()) {
1212 } else {
1214 }
1215 for (auto &s : _serversToDelete) {
1216 delete s;
1217 }
1218 return;
1219 }
1220 }
1221
1222 // generic export using the factory expressions
1223 const auto &dict = exportKeys.find(cl);
1224 if (dict == exportKeys.end()) {
1225 std::cerr << "unable to export class '" << cl->GetName() << "' - no export keys available!\n"
1226 << "there are several possible reasons for this:\n"
1227 << " 1. " << cl->GetName() << " is a custom class that you or some package you are using added.\n"
1228 << " 2. " << cl->GetName()
1229 << " is a ROOT class that nobody ever bothered to write a serialization definition for.\n"
1230 << " 3. something is wrong with your setup, e.g. you might have called "
1231 "RooFit::JSONIO::clearExportKeys() and/or never successfully read a file defining these "
1232 "keys with RooFit::JSONIO::loadExportKeys(filename)\n"
1233 << "either way, please make sure that:\n"
1234 << " 3: you are reading a file with export keys - call RooFit::JSONIO::printExportKeys() to "
1235 "see what is available\n"
1236 << " 2 & 1: you might need to write a serialization definition yourself. check "
1237 "https://root.cern/doc/master/group__roofit__dev__docs__hs3.html to "
1238 "see how to do this!\n";
1239 return;
1240 }
1241
1242 elem["type"] << dict->second.type;
1243
1244 size_t nprox = func.numProxies();
1245
1246 for (size_t i = 0; i < nprox; ++i) {
1247 RooAbsProxy *p = func.getProxy(i);
1248 if (!p)
1249 continue;
1250
1251 // some proxies start with a "!". This is a magic symbol that we don't want to stream
1252 std::string pname(p->name());
1253 if (pname[0] == '!')
1254 pname.erase(0, 1);
1255
1256 auto k = dict->second.proxies.find(pname);
1257 if (k == dict->second.proxies.end()) {
1258 std::cerr << "failed to find key matching proxy '" << pname << "' for type '" << dict->second.type
1259 << "', encountered in '" << func.GetName() << "', skipping" << std::endl;
1260 return;
1261 }
1262
1263 // empty string is interpreted as an instruction to ignore this value
1264 if (k->second.empty())
1265 continue;
1266
1267 if (auto l = dynamic_cast<RooAbsCollection *>(p)) {
1268 fillSeq(elem[k->second], *l);
1269 }
1270 if (auto r = dynamic_cast<RooArgProxy *>(p)) {
1271 if (isLiteralConstVar(*r->absArg())) {
1272 elem[k->second] << static_cast<RooConstVar *>(r->absArg())->getVal();
1273 } else {
1274 elem[k->second] << r->absArg()->GetName();
1275 }
1276 }
1277 }
1278
1279 // export all the servers of a given RooAbsArg
1280 for (RooAbsArg *s : func.servers()) {
1281 if (!s) {
1282 std::cerr << "unable to locate server of " << func.GetName() << std::endl;
1283 continue;
1284 }
1286 }
1287}
1288
1289/**
1290 * @brief Import a function from the JSONNode into the workspace.
1291 *
1292 * This function imports a function from the given JSONNode into the workspace.
1293 * The function's information is read from the JSONNode and added to the workspace.
1294 *
1295 * @param p The JSONNode representing the function to be imported.
1296 * @param importAllDependants A boolean flag indicating whether to import all dependants (servers) of the function.
1297 * @return void
1298 */
1300{
1301 std::string name(RooJSONFactoryWSTool::name(p));
1302
1303 // If this node if marked to be skipped by RooFit, exit
1304 if (hasAttribute(name, "roofit_skip")) {
1305 return;
1306 }
1307
1308 auto const &importers = RooFit::JSONIO::importers();
1310
1311 // some preparations: what type of function are we dealing with here?
1313
1314 // if the RooAbsArg already exists, we don't need to do anything
1315 if (_workspace.arg(name)) {
1316 return;
1317 }
1318 // if the key we found is not a map, it's an error
1319 if (!p.is_map()) {
1320 std::stringstream ss;
1321 ss << "RooJSONFactoryWSTool() function node " + name + " is not a map!";
1323 return;
1324 }
1325 std::string prefix = genPrefix(p, true);
1326 if (!prefix.empty())
1327 name = prefix + name;
1328 if (!p.has_child("type")) {
1329 std::stringstream ss;
1330 ss << "RooJSONFactoryWSTool() no type given for function '" << name << "', skipping." << std::endl;
1332 return;
1333 }
1334
1335 std::string functype(p["type"].val());
1336
1337 // import all dependents if importing a workspace, not for creating new objects
1338 if (!importAllDependants) {
1339 this->importDependants(p);
1340 }
1341
1342 // check for specific implementations
1343 auto it = importers.find(functype);
1344 bool ok = false;
1345 if (it != importers.end()) {
1346 for (auto &imp : it->second) {
1347 try {
1348 ok = imp->importArg(this, p);
1349 } catch (const std::exception &e) {
1350 std::stringstream ss;
1351 const auto *ptr = imp.get();
1352 ss << "RooJSONFactoryWSTool() failed. The importer " << typeid(*ptr).name()
1353 << " emitted and error: " << e.what() << std::endl;
1355 }
1356 if (ok)
1357 break;
1358 }
1359 }
1360 if (!ok) { // generic import using the factory expressions
1361 auto expr = factoryExpressions.find(functype);
1362 if (expr != factoryExpressions.end()) {
1363 std::string expression = ::generate(expr->second, p, this);
1364 if (!_workspace.factory(expression)) {
1365 std::stringstream ss;
1366 ss << "RooJSONFactoryWSTool() failed to create " << expr->second.tclass->GetName() << " '" << name
1367 << "', skipping. expression was\n"
1368 << expression << std::endl;
1370 }
1371 } else {
1372 std::stringstream ss;
1373 ss << "RooJSONFactoryWSTool() no handling for type '" << functype << "' implemented, skipping."
1374 << "\n"
1375 << "there are several possible reasons for this:\n"
1376 << " 1. " << functype << " is a custom type that is not available in RooFit.\n"
1377 << " 2. " << functype
1378 << " is a ROOT class that nobody ever bothered to write a deserialization definition for.\n"
1379 << " 3. something is wrong with your setup, e.g. you might have called "
1380 "RooFit::JSONIO::clearFactoryExpressions() and/or never successfully read a file defining "
1381 "these expressions with RooFit::JSONIO::loadFactoryExpressions(filename)\n"
1382 << "either way, please make sure that:\n"
1383 << " 3: you are reading a file with factory expressions - call "
1384 "RooFit::JSONIO::printFactoryExpressions() "
1385 "to see what is available\n"
1386 << " 2 & 1: you might need to write a deserialization definition yourself. check "
1387 "https://root.cern/doc/master/group__roofit__dev__docs__hs3.html to see "
1388 "how to do this!"
1389 << std::endl;
1391 return;
1392 }
1393 }
1395 if (!func) {
1396 std::stringstream err;
1397 err << "something went wrong importing function '" << name << "'.";
1398 RooJSONFactoryWSTool::error(err.str());
1399 }
1400}
1401
1402/**
1403 * @brief Import a function from a JSON string into the workspace.
1404 *
1405 * This function imports a function from the provided JSON string into the workspace.
1406 * The function's information is read from the JSON string and added to the workspace.
1407 *
1408 * @param jsonString The JSON string containing the function information.
1409 * @param importAllDependants A boolean flag indicating whether to import all dependants (servers) of the function.
1410 * @return void
1411 */
1413{
1414 this->importFunction((JSONTree::create(jsonString))->rootnode(), importAllDependants);
1415}
1416
1417/**
1418 * @brief Export histogram data to a JSONNode.
1419 *
1420 * This function exports histogram data, represented by the provided variables and contents, to a JSONNode.
1421 * The histogram's axes information and bin contents are added as key-value pairs to the JSONNode.
1422 *
1423 * @param vars The RooArgSet representing the variables associated with the histogram.
1424 * @param n The number of bins in the histogram.
1425 * @param contents A pointer to the array containing the bin contents of the histogram.
1426 * @param output The JSONNode to which the histogram data will be exported.
1427 * @return void
1428 */
1429void RooJSONFactoryWSTool::exportHisto(RooArgSet const &vars, std::size_t n, double const *contents, JSONNode &output)
1430{
1431 auto &observablesNode = output["axes"].set_seq();
1432 // axes have to be ordered to get consistent bin indices
1433 for (auto *var : static_range_cast<RooRealVar *>(vars)) {
1434 std::string name = var->GetName();
1436 JSONNode &obsNode = observablesNode.append_child().set_map();
1437 obsNode["name"] << name;
1438 if (var->getBinning().isUniform()) {
1439 obsNode["min"] << var->getMin();
1440 obsNode["max"] << var->getMax();
1441 obsNode["nbins"] << var->getBins();
1442 } else {
1443 auto &edges = obsNode["edges"];
1444 edges.set_seq();
1445 double val = var->getBinning().binLow(0);
1446 edges.append_child() << val;
1447 for (int i = 0; i < var->getBinning().numBins(); ++i) {
1448 val = var->getBinning().binHigh(i);
1449 edges.append_child() << val;
1450 }
1451 }
1452 }
1453
1454 return exportArray(n, contents, output["contents"]);
1455}
1456
1457/**
1458 * @brief Export an array of doubles to a JSONNode.
1459 *
1460 * This function exports an array of doubles, represented by the provided size and contents,
1461 * to a JSONNode. The array elements are added to the JSONNode as a sequence of values.
1462 *
1463 * @param n The size of the array.
1464 * @param contents A pointer to the array containing the double values.
1465 * @param output The JSONNode to which the array will be exported.
1466 * @return void
1467 */
1468void RooJSONFactoryWSTool::exportArray(std::size_t n, double const *contents, JSONNode &output)
1469{
1470 output.set_seq();
1471 for (std::size_t i = 0; i < n; ++i) {
1472 double w = contents[i];
1473 // To make sure there are no unnecessary floating points in the JSON
1474 if (int(w) == w) {
1475 output.append_child() << int(w);
1476 } else {
1477 output.append_child() << w;
1478 }
1479 }
1480}
1481
1482/**
1483 * @brief Export a RooAbsCategory object to a JSONNode.
1484 *
1485 * This function exports a RooAbsCategory object, represented by the provided categories and indices,
1486 * to a JSONNode. The category labels and corresponding indices are added to the JSONNode as key-value pairs.
1487 *
1488 * @param cat The RooAbsCategory object to be exported.
1489 * @param node The JSONNode to which the category data will be exported.
1490 * @return void
1491 */
1493{
1494 auto &labels = node["labels"].set_seq();
1495 auto &indices = node["indices"].set_seq();
1496
1497 for (auto const &item : cat) {
1498 std::string label;
1499 if (std::isalpha(item.first[0])) {
1501 if (label != item.first) {
1502 oocoutW(nullptr, IO) << "RooFitHS3: changed '" << item.first << "' to '" << label
1503 << "' to become a valid name";
1504 }
1505 } else {
1506 RooJSONFactoryWSTool::error("refusing to change first character of string '" + item.first +
1507 "' to make a valid name!");
1508 label = item.first;
1509 }
1510 labels.append_child() << label;
1511 indices.append_child() << item.second;
1512 }
1513}
1514
1515/**
1516 * @brief Export combined data from the workspace to a custom struct.
1517 *
1518 * This function exports combined data from the workspace, represented by the provided RooAbsData object,
1519 * to a CombinedData struct. The struct contains information such as variables, categories,
1520 * and bin contents of the combined data.
1521 *
1522 * @param data The RooAbsData object representing the combined data to be exported.
1523 * @return CombinedData A custom struct containing the exported combined data.
1524 */
1526{
1527 // find category observables
1528 RooAbsCategory *cat = nullptr;
1529 for (RooAbsArg *obs : *data.get()) {
1530 if (dynamic_cast<RooAbsCategory *>(obs)) {
1531 if (cat) {
1532 RooJSONFactoryWSTool::error("dataset '" + std::string(data.GetName()) +
1533 " has several category observables!");
1534 }
1535 cat = static_cast<RooAbsCategory *>(obs);
1536 }
1537 }
1538
1539 // prepare return value
1541
1542 if (!cat)
1543 return datamap;
1544 // this is a combined dataset
1545
1546 datamap.name = data.GetName();
1547
1548 // Write information necessary to reconstruct the combined dataset upon import
1549 auto &child = getRooFitInternal(*_rootnodeOutput, "combined_datasets").set_map()[data.GetName()].set_map();
1550 child["index_cat"] << cat->GetName();
1551 exportCategory(*cat, child);
1552
1553 // Find a RooSimultaneous model that would fit to this dataset
1554 RooSimultaneous const *simPdf = nullptr;
1555 auto *combinedPdfInfoNode = findRooFitInternal(*_rootnodeOutput, "combined_distributions");
1556 if (combinedPdfInfoNode) {
1557 for (auto &info : combinedPdfInfoNode->children()) {
1558 if (info["index_cat"].val() == cat->GetName()) {
1559 simPdf = static_cast<RooSimultaneous const *>(_workspace.pdf(info.key()));
1560 }
1561 }
1562 }
1563
1564 // If there is an associated simultaneous pdf for the index category, we
1565 // use the RooAbsData::split() overload that takes the RooSimultaneous.
1566 // Like this, the observables that are not relevant for a given channel
1567 // are automatically split from the component datasets.
1568 std::vector<std::unique_ptr<RooAbsData>> dataList{simPdf ? data.split(*simPdf, true) : data.split(*cat, true)};
1569
1570 for (std::unique_ptr<RooAbsData> const &absData : dataList) {
1571 std::string catName(absData->GetName());
1572 std::string dataName;
1573 if (std::isalpha(catName[0])) {
1575 if (dataName != catName) {
1576 oocoutW(nullptr, IO) << "RooFitHS3: changed '" << catName << "' to '" << dataName
1577 << "' to become a valid name";
1578 }
1579 } else {
1580 RooJSONFactoryWSTool::error("refusing to change first character of string '" + catName +
1581 "' to make a valid name!");
1582 dataName = catName;
1583 }
1584 absData->SetName((std::string(data.GetName()) + "_" + dataName).c_str());
1585 datamap.components[catName] = absData->GetName();
1586 this->exportData(*absData);
1587 }
1588 return datamap;
1589}
1590
1591/**
1592 * @brief Export data from the workspace to a JSONNode.
1593 *
1594 * This function exports data represented by the provided RooAbsData object,
1595 * to a JSONNode. The data's information is added as key-value pairs to the JSONNode.
1596 *
1597 * @param data The RooAbsData object representing the data to be exported.
1598 * @return void
1599 */
1601{
1602 // find category observables
1603
1604 RooAbsCategory *cat = nullptr;
1605 for (RooAbsArg *obs : *data.get()) {
1606 if (dynamic_cast<RooAbsCategory *>(obs)) {
1607 if (cat) {
1608 RooJSONFactoryWSTool::error("dataset '" + std::string(data.GetName()) +
1609 " has several category observables!");
1610 }
1611 cat = static_cast<RooAbsCategory *>(obs);
1612 }
1613 }
1614
1615 if (cat)
1616 return;
1617
1618 JSONNode &output = appendNamedChild((*_rootnodeOutput)["data"], data.GetName());
1619
1620 // This works around a problem in RooStats/HistFactory that was only fixed
1621 // in ROOT 6.30: until then, the weight variable of the observed dataset,
1622 // called "weightVar", was added to the observables. Therefore, it also got
1623 // added to the Asimov dataset. But the Asimov has its own weight variable,
1624 // called "binWeightAsimov", making "weightVar" an actual observable in the
1625 // Asimov data. But this is only by accident and should be removed.
1626 RooArgSet variables = *data.get();
1627 if (auto weightVar = variables.find("weightVar")) {
1628 variables.remove(*weightVar);
1629 }
1630
1631 // this is a regular binned dataset
1632 if (auto dh = dynamic_cast<RooDataHist const *>(&data)) {
1633 output["type"] << "binned";
1634 for (auto *var : static_range_cast<RooRealVar *>(variables)) {
1635 _domains->readVariable(*var);
1636 }
1637 return exportHisto(variables, dh->numEntries(), dh->weightArray(), output);
1638 }
1639
1640 // Check if this actually represents a binned dataset, and then import it
1641 // like a RooDataHist. This happens frequently when people create combined
1642 // RooDataSets from binned data to fit HistFactory models. In this case, it
1643 // doesn't make sense to export them like an unbinned dataset, because the
1644 // coordinates are redundant information with the binning. We only do this
1645 // for 1D data for now.
1646 if (data.isWeighted() && variables.size() == 1) {
1647 bool isBinnedData = false;
1648 auto &x = static_cast<RooRealVar const &>(*variables[0]);
1649 std::vector<double> contents;
1650 int i = 0;
1651 for (; i < data.numEntries(); ++i) {
1652 data.get(i);
1653 if (x.getBin() != i)
1654 break;
1655 contents.push_back(data.weight());
1656 }
1657 if (i == x.getBins())
1658 isBinnedData = true;
1659 if (isBinnedData) {
1660 output["type"] << "binned";
1661 for (auto *var : static_range_cast<RooRealVar *>(variables)) {
1662 _domains->readVariable(*var);
1663 }
1664 return exportHisto(variables, data.numEntries(), contents.data(), output);
1665 }
1666 }
1667
1668 // this really is an unbinned dataset
1669 output["type"] << "unbinned";
1670 exportVariables(variables, output["axes"], false, true);
1671 auto &coords = output["entries"].set_seq();
1672 std::vector<double> weightVals;
1673 bool hasNonUnityWeights = false;
1674 for (int i = 0; i < data.numEntries(); ++i) {
1675 data.get(i);
1676 coords.append_child().fill_seq(variables, [](auto x) { return static_cast<RooRealVar *>(x)->getVal(); });
1677 std::string datasetName = data.GetName();
1678 if (data.isWeighted()) {
1679 weightVals.push_back(data.weight());
1680 if (data.weight() != 1.)
1681 hasNonUnityWeights = true;
1682 }
1683 }
1684 if (data.isWeighted() && hasNonUnityWeights) {
1685 output["weights"].fill_seq(weightVals);
1686 }
1687}
1688
1689/**
1690 * @brief Read axes from the JSONNode and create a RooArgSet representing them.
1691 *
1692 * This function reads axes information from the given JSONNode and
1693 * creates a RooArgSet with variables representing these axes.
1694 *
1695 * @param topNode The JSONNode containing the axes information to be read.
1696 * @return RooArgSet A RooArgSet containing the variables created from the JSONNode.
1697 */
1699{
1700 RooArgSet vars;
1701
1702 for (JSONNode const &node : topNode["axes"].children()) {
1703 if (node.has_child("edges")) {
1704 std::vector<double> edges;
1705 for (auto const &bound : node["edges"].children()) {
1706 edges.push_back(bound.val_double());
1707 }
1708 auto obs = std::make_unique<RooRealVar>(node["name"].val().c_str(), node["name"].val().c_str(), edges[0],
1709 edges[edges.size() - 1]);
1710 RooBinning bins(obs->getMin(), obs->getMax());
1711 for (auto b : edges) {
1712 bins.addBoundary(b);
1713 }
1714 obs->setBinning(bins);
1715 vars.addOwned(std::move(obs));
1716 } else {
1717 auto obs = std::make_unique<RooRealVar>(node["name"].val().c_str(), node["name"].val().c_str(),
1718 node["min"].val_double(), node["max"].val_double());
1719 obs->setBins(node["nbins"].val_int());
1720 vars.addOwned(std::move(obs));
1721 }
1722 }
1723
1724 return vars;
1725}
1726
1727/**
1728 * @brief Read binned data from the JSONNode and create a RooDataHist object.
1729 *
1730 * This function reads binned data from the given JSONNode and creates a RooDataHist object.
1731 * The binned data is associated with the specified name and variables (RooArgSet) in the workspace.
1732 *
1733 * @param n The JSONNode representing the binned data to be read.
1734 * @param name The name to be associated with the created RooDataHist object.
1735 * @param vars The RooArgSet representing the variables associated with the binned data.
1736 * @return std::unique_ptr<RooDataHist> A unique pointer to the created RooDataHist object.
1737 */
1738std::unique_ptr<RooDataHist>
1739RooJSONFactoryWSTool::readBinnedData(const JSONNode &n, const std::string &name, RooArgSet const &vars)
1740{
1741 if (!n.has_child("contents"))
1742 RooJSONFactoryWSTool::error("no contents given");
1743
1744 JSONNode const &contents = n["contents"];
1745
1746 if (!contents.is_seq())
1747 RooJSONFactoryWSTool::error("contents are not in list form");
1748
1749 JSONNode const *errors = nullptr;
1750 if (n.has_child("errors")) {
1751 errors = &n["errors"];
1752 if (!errors->is_seq())
1753 RooJSONFactoryWSTool::error("errors are not in list form");
1754 }
1755
1756 auto bins = generateBinIndices(vars);
1757 if (contents.num_children() != bins.size()) {
1758 std::stringstream errMsg;
1759 errMsg << "inconsistent bin numbers: contents=" << contents.num_children() << ", bins=" << bins.size();
1761 }
1762 auto dh = std::make_unique<RooDataHist>(name, name, vars);
1763 std::vector<double> contentVals;
1764 contentVals.reserve(contents.num_children());
1765 for (auto const &cont : contents.children()) {
1766 contentVals.push_back(cont.val_double());
1767 }
1768 std::vector<double> errorVals;
1769 if (errors) {
1770 errorVals.reserve(errors->num_children());
1771 for (auto const &err : errors->children()) {
1772 errorVals.push_back(err.val_double());
1773 }
1774 }
1775 for (size_t ibin = 0; ibin < bins.size(); ++ibin) {
1776 const double err = errors ? errorVals[ibin] : -1;
1777 dh->set(ibin, contentVals[ibin], err);
1778 }
1779 return dh;
1780}
1781
1782/**
1783 * @brief Import a variable from the JSONNode into the workspace.
1784 *
1785 * This function imports a variable from the given JSONNode into the workspace.
1786 * The variable's information is read from the JSONNode and added to the workspace.
1787 *
1788 * @param p The JSONNode representing the variable to be imported.
1789 * @return void
1790 */
1792{
1793 // import a RooRealVar object
1794 std::string name(RooJSONFactoryWSTool::name(p));
1796
1797 if (_workspace.var(name))
1798 return;
1799 if (!p.is_map()) {
1800 std::stringstream ss;
1801 ss << "RooJSONFactoryWSTool() node '" << name << "' is not a map, skipping.";
1802 oocoutE(nullptr, InputArguments) << ss.str() << std::endl;
1803 return;
1804 }
1805 if (_attributesNode) {
1806 if (auto *attrNode = _attributesNode->find(name)) {
1807 // We should not create RooRealVar objects for RooConstVars!
1808 if (attrNode->has_child("is_const_var") && (*attrNode)["is_const_var"].val_int() == 1) {
1809 wsEmplace<RooConstVar>(name, p["value"].val_double());
1810 return;
1811 }
1812 }
1813 }
1815}
1816
1817/**
1818 * @brief Import all dependants (servers) of a node into the workspace.
1819 *
1820 * This function imports all the dependants (servers) of the given JSONNode into the workspace.
1821 * The dependants' information is read from the JSONNode and added to the workspace.
1822 *
1823 * @param n The JSONNode representing the node whose dependants are to be imported.
1824 * @return void
1825 */
1827{
1828 // import all the dependants of an object
1829 if (JSONNode const *varsNode = getVariablesNode(n)) {
1830 for (const auto &p : varsNode->children()) {
1832 }
1833 }
1834 if (auto seq = n.find("functions")) {
1835 for (const auto &p : seq->children()) {
1836 this->importFunction(p, true);
1837 }
1838 }
1839 if (auto seq = n.find("distributions")) {
1840 for (const auto &p : seq->children()) {
1841 this->importFunction(p, true);
1842 }
1843 }
1844}
1845
1847 const std::vector<CombinedData> &combDataSets,
1848 const std::vector<RooAbsData *> &singleDataSets)
1849{
1850 auto pdf = mc.GetPdf();
1851 auto simpdf = dynamic_cast<RooSimultaneous const *>(pdf);
1852 if (simpdf) {
1853 for (std::size_t i = 0; i < std::max(combDataSets.size(), std::size_t(1)); ++i) {
1854 const bool hasdata = i < combDataSets.size();
1855 if (hasdata && !matches(combDataSets.at(i), simpdf))
1856 continue;
1857
1858 std::string analysisName(simpdf->GetName());
1859 if (hasdata)
1860 analysisName += "_" + combDataSets[i].name;
1861
1862 exportSingleModelConfig(rootnode, mc, analysisName, hasdata ? &combDataSets[i].components : nullptr);
1863 }
1864 } else {
1865 RooArgSet observables(*mc.GetObservables());
1866 int founddata = 0;
1867 for (auto *data : singleDataSets) {
1868 if (observables.equals(*(data->get()))) {
1869 std::map<std::string, std::string> mapping;
1870 mapping[pdf->GetName()] = data->GetName();
1871 exportSingleModelConfig(rootnode, mc, std::string(pdf->GetName()) + "_" + data->GetName(), &mapping);
1872 ++founddata;
1873 }
1874 }
1875 if (founddata == 0) {
1876 exportSingleModelConfig(rootnode, mc, pdf->GetName(), nullptr);
1877 }
1878 }
1879}
1880
1882 std::string const &analysisName,
1883 std::map<std::string, std::string> const *dataComponents)
1884{
1885 auto pdf = mc.GetPdf();
1886
1887 JSONNode &analysisNode = appendNamedChild(rootnode["analyses"], analysisName);
1888
1889 auto &domains = analysisNode["domains"].set_seq();
1890
1891 analysisNode["likelihood"] << analysisName;
1892
1893 auto &nllNode = appendNamedChild(rootnode["likelihoods"], analysisName);
1894 nllNode["distributions"].set_seq();
1895 nllNode["data"].set_seq();
1896
1897 if (dataComponents) {
1898 auto simPdf = dynamic_cast<RooSimultaneous const *>(pdf);
1899 if (simPdf) {
1900 for (auto const &item : simPdf->indexCat()) {
1901 const auto &dataComp = dataComponents->find(item.first);
1902 nllNode["distributions"].append_child() << simPdf->getPdf(item.first)->GetName();
1903 nllNode["data"].append_child() << dataComp->second;
1904 }
1905 } else {
1906 for (auto it : *dataComponents) {
1907 nllNode["distributions"].append_child() << it.first;
1908 nllNode["data"].append_child() << it.second;
1909 }
1910 }
1911 } else {
1912 nllNode["distributions"].append_child() << pdf->GetName();
1913 nllNode["data"].append_child() << 0;
1914 }
1915
1916 if (mc.GetExternalConstraints()) {
1917 auto &extConstrNode = nllNode["aux_distributions"];
1918 extConstrNode.set_seq();
1919 for (const auto &constr : *mc.GetExternalConstraints()) {
1920 extConstrNode.append_child() << constr->GetName();
1921 }
1922 }
1923
1924 auto writeList = [&](const char *name, RooArgSet const *args) {
1925 if (!args || !args->size())
1926 return;
1927
1928 std::vector<std::string> names;
1929 names.reserve(args->size());
1930 for (RooAbsArg const *arg : *args)
1931 names.push_back(arg->GetName());
1932 std::sort(names.begin(), names.end());
1933 analysisNode[name].fill_seq(names);
1934 };
1935
1936 writeList("parameters_of_interest", mc.GetParametersOfInterest());
1937
1938 auto &domainsNode = rootnode["domains"];
1939
1940 auto writeProductDomain = [&](const char *suffix, RooArgSet const *args) {
1941 if (!args || args->empty())
1942 return;
1943 const std::string domainName = analysisName + suffix;
1944 domains.append_child() << domainName;
1946 for (auto *var : static_range_cast<const RooRealVar *>(*args)) {
1947 domain.readVariable(*var);
1948 }
1950 };
1951
1952 writeProductDomain("_nuisance_parameters", mc.GetNuisanceParameters());
1953 writeProductDomain("_global_observables", mc.GetGlobalObservables());
1954 writeProductDomain("_parameters_of_interest", mc.GetParametersOfInterest());
1955
1956 auto &modelConfigAux = getRooFitInternal(rootnode, "ModelConfigs", analysisName);
1957 modelConfigAux.set_map();
1958 modelConfigAux["pdfName"] << pdf->GetName();
1959 modelConfigAux["mcName"] << mc.GetName();
1960}
1961
1962/**
1963 * @brief Export all objects in the workspace to a JSONNode.
1964 *
1965 * This function exports all the objects in the workspace to the provided JSONNode.
1966 * The objects' information is added as key-value pairs to the JSONNode.
1967 *
1968 * @param n The JSONNode to which the objects will be exported.
1969 * @return void
1970 */
1972{
1973 _domains = std::make_unique<RooFit::JSONIO::Detail::Domains>();
1975 _rootnodeOutput = &n;
1976
1977 // export all toplevel pdfs
1978 std::vector<RooAbsPdf *> allpdfs;
1979 for (auto &arg : _workspace.allPdfs()) {
1980 if (!arg->hasClients()) {
1981 if (auto *pdf = dynamic_cast<RooAbsPdf *>(arg)) {
1982 allpdfs.push_back(pdf);
1983 }
1984 }
1985 }
1987 std::set<std::string> exportedObjectNames;
1989
1990 // export all toplevel functions
1991 std::vector<RooAbsReal *> allfuncs;
1992 for (auto &arg : _workspace.allFunctions()) {
1993 if (!arg->hasClients()) {
1994 if (auto *func = dynamic_cast<RooAbsReal *>(arg)) {
1995 allfuncs.push_back(func);
1996 }
1997 }
1998 }
2001
2002 // export attributes of all objects
2003 for (RooAbsArg *arg : _workspace.components()) {
2004 exportAttributes(arg, n);
2005 }
2006
2007 // collect all datasets
2008 std::vector<RooAbsData *> alldata;
2009 for (auto &d : _workspace.allData()) {
2010 alldata.push_back(d);
2011 }
2013 // first, take care of combined datasets
2014 std::vector<RooAbsData *> singleData;
2015 std::vector<RooJSONFactoryWSTool::CombinedData> combData;
2016 for (auto &d : alldata) {
2017 auto data = this->exportCombinedData(*d);
2018 if (!data.components.empty())
2019 combData.push_back(data);
2020 else
2021 singleData.push_back(d);
2022 }
2023 // next, take care datasets
2024 for (auto &d : alldata) {
2025 this->exportData(*d);
2026 }
2027
2028 // export all ModelConfig objects and attached Pdfs
2029 for (TObject *obj : _workspace.allGenericObjects()) {
2030 if (auto mc = dynamic_cast<RooFit::ModelConfig *>(obj)) {
2032 }
2033 }
2034
2036
2039 // We only want to add the variables that actually got exported and skip
2040 // the ones that the pdfs encoded implicitly (like in the case of
2041 // HistFactory).
2042 for (RooAbsArg *arg : *snsh) {
2043 bool do_export = false;
2044 for (const auto &pdf : allpdfs) {
2045 if (pdf->dependsOn(*arg)) {
2046 do_export = true;
2047 }
2048 }
2049 if (do_export) {
2050 RooJSONFactoryWSTool::testValidName(arg->GetName(), true);
2051 snapshotSorted.add(*arg);
2052 }
2053 }
2054 snapshotSorted.sort();
2055 std::string name(snsh->GetName());
2056 if (name != "default_values") {
2057 this->exportVariables(snapshotSorted, appendNamedChild(n["parameter_points"], name)["parameters"], true,
2058 false);
2059 }
2060 }
2061 _varsNode = nullptr;
2062 _domains->writeJSON(n["domains"]);
2063 _domains.reset();
2064 _rootnodeOutput = nullptr;
2065}
2066
2067/**
2068 * @brief Import the workspace from a JSON string.
2069 *
2070 * @param s The JSON string containing the workspace data.
2071 * @return bool Returns true on successful import, false otherwise.
2072 */
2074{
2075 std::stringstream ss(s);
2076 return importJSON(ss);
2077}
2078
2079/**
2080 * @brief Import the workspace from a YML string.
2081 *
2082 * @param s The YML string containing the workspace data.
2083 * @return bool Returns true on successful import, false otherwise.
2084 */
2086{
2087 std::stringstream ss(s);
2088 return importYML(ss);
2089}
2090
2091/**
2092 * @brief Export the workspace to a JSON string.
2093 *
2094 * @return std::string The JSON string representing the exported workspace.
2095 */
2097{
2098 std::stringstream ss;
2099 exportJSON(ss);
2100 return ss.str();
2101}
2102
2103/**
2104 * @brief Export the workspace to a YML string.
2105 *
2106 * @return std::string The YML string representing the exported workspace.
2107 */
2109{
2110 std::stringstream ss;
2111 exportYML(ss);
2112 return ss.str();
2113}
2114
2115/**
2116 * @brief Create a new JSON tree with version information.
2117 *
2118 * @return std::unique_ptr<JSONTree> A unique pointer to the created JSON tree.
2119 */
2121{
2122 std::unique_ptr<JSONTree> tree = JSONTree::create();
2123 JSONNode &n = tree->rootnode();
2124 n.set_map();
2125 auto &metadata = n["metadata"].set_map();
2126
2127 // add the mandatory hs3 version number
2128 metadata["hs3_version"] << hs3VersionTag;
2129
2130 // Add information about the ROOT version that was used to generate this file
2131 auto &rootInfo = appendNamedChild(metadata["packages"], "ROOT");
2132 std::string versionName = gROOT->GetVersion();
2133 // We want to consistently use dots such that the version name can be easily
2134 // digested automatically.
2135 std::replace(versionName.begin(), versionName.end(), '/', '.');
2136 rootInfo["version"] << versionName;
2137
2138 return tree;
2139}
2140
2141/**
2142 * @brief Export the workspace to JSON format and write to the output stream.
2143 *
2144 * @param os The output stream to write the JSON data to.
2145 * @return bool Returns true on successful export, false otherwise.
2146 */
2148{
2149 std::unique_ptr<JSONTree> tree = createNewJSONTree();
2150 JSONNode &n = tree->rootnode();
2151 this->exportAllObjects(n);
2152 n.writeJSON(os);
2153 return true;
2154}
2155
2156/**
2157 * @brief Export the workspace to JSON format and write to the specified file.
2158 *
2159 * @param filename The name of the JSON file to create and write the data to.
2160 * @return bool Returns true on successful export, false otherwise.
2161 */
2163{
2164 std::ofstream out(filename.c_str());
2165 if (!out.is_open()) {
2166 std::stringstream ss;
2167 ss << "RooJSONFactoryWSTool() invalid output file '" << filename << "'." << std::endl;
2169 return false;
2170 }
2171 return this->exportJSON(out);
2172}
2173
2174/**
2175 * @brief Export the workspace to YML format and write to the output stream.
2176 *
2177 * @param os The output stream to write the YML data to.
2178 * @return bool Returns true on successful export, false otherwise.
2179 */
2181{
2182 std::unique_ptr<JSONTree> tree = createNewJSONTree();
2183 JSONNode &n = tree->rootnode();
2184 this->exportAllObjects(n);
2185 n.writeYML(os);
2186 return true;
2187}
2188
2189/**
2190 * @brief Export the workspace to YML format and write to the specified file.
2191 *
2192 * @param filename The name of the YML file to create and write the data to.
2193 * @return bool Returns true on successful export, false otherwise.
2194 */
2196{
2197 std::ofstream out(filename.c_str());
2198 if (!out.is_open()) {
2199 std::stringstream ss;
2200 ss << "RooJSONFactoryWSTool() invalid output file '" << filename << "'." << std::endl;
2202 return false;
2203 }
2204 return this->exportYML(out);
2205}
2206
2207bool RooJSONFactoryWSTool::hasAttribute(const std::string &obj, const std::string &attrib)
2208{
2209 if (!_attributesNode)
2210 return false;
2211 if (auto attrNode = _attributesNode->find(obj)) {
2212 if (auto seq = attrNode->find("tags")) {
2213 for (auto &a : seq->children()) {
2214 if (a.val() == attrib)
2215 return true;
2216 }
2217 }
2218 }
2219 return false;
2220}
2221void RooJSONFactoryWSTool::setAttribute(const std::string &obj, const std::string &attrib)
2222{
2223 auto node = &RooJSONFactoryWSTool::getRooFitInternal(*_rootnodeOutput, "attributes").set_map()[obj].set_map();
2224 auto &tags = (*node)["tags"];
2225 tags.set_seq();
2226 tags.append_child() << attrib;
2227}
2228
2229std::string RooJSONFactoryWSTool::getStringAttribute(const std::string &obj, const std::string &attrib)
2230{
2231 if (!_attributesNode)
2232 return "";
2233 if (auto attrNode = _attributesNode->find(obj)) {
2234 if (auto dict = attrNode->find("dict")) {
2235 if (auto *a = dict->find(attrib)) {
2236 return a->val();
2237 }
2238 }
2239 }
2240 return "";
2241}
2242void RooJSONFactoryWSTool::setStringAttribute(const std::string &obj, const std::string &attrib,
2243 const std::string &value)
2244{
2245 auto node = &RooJSONFactoryWSTool::getRooFitInternal(*_rootnodeOutput, "attributes").set_map()[obj].set_map();
2246 auto &dict = (*node)["dict"];
2247 dict.set_map();
2248 dict[attrib] << value;
2249}
2250
2251/**
2252 * @brief Imports all nodes of the JSON data and adds them to the workspace.
2253 *
2254 * @param n The JSONNode representing the root node of the JSON data.
2255 * @return void
2256 */
2258{
2259 // Per HS3 standard, the hs3_version in the metadata is required. So we
2260 // error out if it is missing. TODO: now we are only checking if the
2261 // hs3_version tag exists, but in the future when the HS3 specification
2262 // versions are actually frozen, we should also check if the hs3_version is
2263 // one that RooFit can actually read.
2264 auto metadata = n.find("metadata");
2265 if (!metadata || !metadata->find("hs3_version")) {
2266 std::stringstream ss;
2267 ss << "The HS3 version is missing in the JSON!\n"
2268 << "Please include the HS3 version in the metadata field, e.g.:\n"
2269 << " \"metadata\" :\n"
2270 << " {\n"
2271 << " \"hs3_version\" : \"" << hs3VersionTag << "\"\n"
2272 << " }";
2273 error(ss.str());
2274 }
2275
2276 _domains = std::make_unique<RooFit::JSONIO::Detail::Domains>();
2277 if (auto domains = n.find("domains")) {
2278 _domains->readJSON(*domains);
2279 }
2280 _domains->populate(_workspace);
2281
2282 _rootnodeInput = &n;
2283
2285
2286 // Build name-keyed indices over the "functions" and "distributions"
2287 // sequences. Without these, every cross-reference resolved during import
2288 // (e.g. dependencies of a PiecewiseInterpolation, or factory-expression
2289 // arguments) triggers a linear scan over all sibling nodes via
2290 // findNamedChild(), which becomes O(N^2) on workspaces with thousands of
2291 // entries. Populating the maps up-front turns each lookup into O(1).
2292 _functionsByName.clear();
2293 _distributionsByName.clear();
2294 if (auto seq = n.find("functions")) {
2295 if (seq->is_seq()) {
2296 _functionsByName.reserve(seq->num_children());
2297 for (const auto &p : seq->children()) {
2299 }
2300 }
2301 }
2302 if (auto seq = n.find("distributions")) {
2303 if (seq->is_seq()) {
2304 _distributionsByName.reserve(seq->num_children());
2305 for (const auto &p : seq->children()) {
2307 }
2308 }
2309 }
2310
2311 this->importDependants(n);
2312
2313 if (auto paramPointsNode = n.find("parameter_points")) {
2314 for (const auto &snsh : paramPointsNode->children()) {
2315 std::string name = RooJSONFactoryWSTool::name(snsh);
2317
2318 RooArgSet vars;
2319 for (const auto &var : snsh["parameters"].children()) {
2322 vars.add(*rrv);
2323 }
2324 }
2326 }
2327 }
2328
2330
2331 // Import attributes
2332 if (_attributesNode) {
2333 for (const auto &elem : _attributesNode->children()) {
2334 if (RooAbsArg *arg = _workspace.arg(elem.key()))
2335 importAttributes(arg, elem);
2336 }
2337 }
2338
2339 _attributesNode = nullptr;
2340
2341 // We delay the import of the data to after combineDatasets(), because it
2342 // might be that some datasets are merged to combined datasets there. In
2343 // that case, we will remove the components from the "datasets" vector so they
2344 // don't get imported.
2345 std::vector<std::unique_ptr<RooAbsData>> datasets;
2346 if (auto dataNode = n.find("data")) {
2347 for (const auto &p : dataNode->children()) {
2348 datasets.push_back(loadData(p, _workspace));
2349 }
2350 }
2351
2352 // Now, read in analyses and likelihoods if there are any
2353
2354 if (auto analysesNode = n.find("analyses")) {
2355 for (JSONNode const &analysisNode : analysesNode->children()) {
2356 importAnalysis(*_rootnodeInput, analysisNode, n["likelihoods"], n["domains"], _workspace, datasets);
2357 }
2358 }
2359
2360 combineDatasets(*_rootnodeInput, datasets);
2361
2362 for (auto const &d : datasets) {
2363 if (d) {
2365 for (auto const &obs : *d->get()) {
2366 if (auto *rrv = dynamic_cast<RooRealVar *>(obs)) {
2367 _workspace.var(rrv->GetName())->setBinning(rrv->getBinning());
2368 }
2369 }
2370 }
2371 }
2372
2373 _rootnodeInput = nullptr;
2374 _domains.reset();
2375 _functionsByName.clear();
2376 _distributionsByName.clear();
2377}
2378
2379/**
2380 * @brief Imports a JSON file from the given input stream to the workspace.
2381 *
2382 * @param is The input stream containing the JSON data.
2383 * @return bool Returns true on successful import, false otherwise.
2384 */
2386{
2387 // import a JSON file to the workspace
2388 std::unique_ptr<JSONTree> tree = JSONTree::create(is);
2389 JSONNode const &rootnode = tree->rootnode();
2390 this->importAllNodes(rootnode);
2391 if (this->workspace()->getSnapshot("default_values")) {
2392 this->workspace()->loadSnapshot("default_values");
2393 }
2394 importParameterStepWidths(*this->workspace(), rootnode);
2395 return true;
2396}
2397
2398/**
2399 * @brief Imports a JSON file from the given filename to the workspace.
2400 *
2401 * @param filename The name of the JSON file to import.
2402 * @return bool Returns true on successful import, false otherwise.
2403 */
2405{
2406 // import a JSON file to the workspace
2407 std::ifstream infile(filename.c_str());
2408 if (!infile.is_open()) {
2409 std::stringstream ss;
2410 ss << "RooJSONFactoryWSTool() invalid input file '" << filename << "'." << std::endl;
2412 return false;
2413 }
2414 return this->importJSON(infile);
2415}
2416
2417/**
2418 * @brief Imports a YML file from the given input stream to the workspace.
2419 *
2420 * @param is The input stream containing the YML data.
2421 * @return bool Returns true on successful import, false otherwise.
2422 */
2424{
2425 // import a YML file to the workspace
2426 std::unique_ptr<JSONTree> tree = JSONTree::create(is);
2427 JSONNode const &rootnode = tree->rootnode();
2428 this->importAllNodes(rootnode);
2429 importParameterStepWidths(*this->workspace(), rootnode);
2430 return true;
2431}
2432
2433/**
2434 * @brief Imports a YML file from the given filename to the workspace.
2435 *
2436 * @param filename The name of the YML file to import.
2437 * @return bool Returns true on successful import, false otherwise.
2438 */
2440{
2441 // import a YML file to the workspace
2442 std::ifstream infile(filename.c_str());
2443 if (!infile.is_open()) {
2444 std::stringstream ss;
2445 ss << "RooJSONFactoryWSTool() invalid input file '" << filename << "'." << std::endl;
2447 return false;
2448 }
2449 return this->importYML(infile);
2450}
2451
2452void RooJSONFactoryWSTool::importJSONElement(const std::string &name, const std::string &jsonString)
2453{
2454 std::unique_ptr<RooFit::Detail::JSONTree> tree = RooFit::Detail::JSONTree::create(jsonString);
2455 JSONNode &n = tree->rootnode();
2456 n["name"] << name;
2457
2458 bool isVariable = true;
2459 if (n.find("type")) {
2460 isVariable = false;
2461 }
2462
2463 if (isVariable) {
2464 this->importVariableElement(n);
2465 } else {
2466 this->importFunction(n, false);
2467 }
2468}
2469
2471{
2472 std::unique_ptr<RooFit::Detail::JSONTree> tree = varJSONString(elementNode);
2473 JSONNode &n = tree->rootnode();
2474 _domains = std::make_unique<RooFit::JSONIO::Detail::Domains>();
2475 if (auto domains = n.find("domains"))
2476 _domains->readJSON(*domains);
2477
2478 _rootnodeInput = &n;
2480
2482 const auto &p = varsNode->child(0);
2484
2485 auto paramPointsNode = n.find("parameter_points");
2486 const auto &snsh = paramPointsNode->child(0);
2487 std::string name = RooJSONFactoryWSTool::name(snsh);
2488 RooArgSet vars;
2489 const auto &var = snsh["parameters"].child(0);
2492 vars.add(*rrv);
2493 }
2494
2495 // Import attributes
2496 if (_attributesNode) {
2497 for (const auto &elem : _attributesNode->children()) {
2498 if (RooAbsArg *arg = _workspace.arg(elem.key()))
2499 importAttributes(arg, elem);
2500 }
2501 }
2502
2503 _attributesNode = nullptr;
2504 _rootnodeInput = nullptr;
2505 _domains.reset();
2506}
2507
2508/**
2509 * @brief Writes a warning message to the RooFit message service.
2510 *
2511 * @param str The warning message to be logged.
2512 * @return std::ostream& A reference to the output stream.
2513 */
2514std::ostream &RooJSONFactoryWSTool::warning(std::string const &str)
2515{
2516 return RooMsgService::instance().log(nullptr, RooFit::MsgLevel::ERROR, RooFit::IO) << str << std::endl;
2517}
2518
2519/**
2520 * @brief Writes an error message to the RooFit message service and throws a runtime_error.
2521 *
2522 * @param s The error message to be logged and thrown.
2523 * @return void
2524 */
2526{
2527 RooMsgService::instance().log(nullptr, RooFit::MsgLevel::ERROR, RooFit::IO) << s << std::endl;
2528 throw std::runtime_error(s);
2529}
2530
2531/**
2532 * @brief Cleans up names to the HS3 standard
2533 *
2534 * @param str The string to be sanitized.
2535 * @return std::string
2536 */
2537std::string RooJSONFactoryWSTool::sanitizeName(const std::string str)
2538{
2539 std::string result;
2541 for (char c : str) {
2542 switch (c) {
2543 case '[':
2544 case '|':
2545 case ',':
2546 case '(': result += '_'; break;
2547 case ']':
2548 case ')':
2549 // skip these characters entirely
2550 break;
2551 case '.': result += "_dot_"; break;
2552 case '@': result += "at"; break;
2553 case '-': result += "minus"; break;
2554 case '/': result += "_div_"; break;
2555
2556 default: result += c; break;
2557 }
2558 }
2559 return result;
2560 }
2561 return str;
2562}
2563
2565{
2566 // Variables
2567
2569 if (onlyModelConfig) {
2570 for (auto *obj : ws.allGenericObjects()) {
2571 if (auto *mc = dynamic_cast<RooFit::ModelConfig *>(obj)) {
2572 tmpWS.import(*mc->GetPdf(), RooFit::RecycleConflictNodes(true));
2573 }
2574 }
2575
2576 } else {
2577
2578 for (auto *pdf : ws.allPdfs()) {
2579 if (!pdf->hasClients()) {
2580 tmpWS.import(*pdf, RooFit::RecycleConflictNodes(true));
2581 }
2582 }
2583
2584 for (auto *func : ws.allFunctions()) {
2585 if (!func->hasClients()) {
2586 tmpWS.import(*func, RooFit::RecycleConflictNodes(true));
2587 }
2588 }
2589 }
2590
2591 for (auto *data : ws.allData()) {
2592 tmpWS.import(*data);
2593 }
2594
2595 for (auto *obj : ws.allGenericObjects()) {
2596 tmpWS.import(*obj);
2597 }
2598
2599 for (auto *obj : ws.allResolutionModels()) {
2600 tmpWS.import(*obj);
2601 }
2602
2603 for (auto *snsh : ws.getSnapshots()) {
2604 auto *snshSet = dynamic_cast<RooArgSet *>(snsh);
2605 if (snshSet) {
2606 tmpWS.saveSnapshot(snshSet->GetName(), *snshSet, true);
2607 }
2608 }
2609
2610 return tmpWS;
2611}
2612
2613// Sanitize all names in the workspace to be HS3 compliant
2615{
2616 // Variables
2617
2618 RooWorkspace tmpWS = cleanWS(ws, false);
2619
2620 auto sanitizeIfNeeded = [](auto const &list) {
2621 for (auto *obj : list) {
2622 if (!isValidName(obj->GetName())) {
2623 obj->SetName(sanitizeName(obj->GetName()).c_str());
2624 }
2625 }
2626 };
2627 sanitizeIfNeeded(tmpWS.allVars());
2628 sanitizeIfNeeded(tmpWS.allFunctions());
2629 sanitizeIfNeeded(tmpWS.allPdfs());
2630 sanitizeIfNeeded(tmpWS.allResolutionModels());
2631 // Datasets
2632 for (auto *data : tmpWS.allData()) {
2633 // Sanitize dataset name
2634 if (!isValidName(data->GetName())) {
2635 data->SetName(sanitizeName(data->GetName()).c_str());
2636 }
2637 for (auto *obj : *data->get()) {
2638 obj->SetName(sanitizeName(obj->GetName()).c_str());
2639 }
2640 }
2641 for (auto *data : tmpWS.allEmbeddedData()) {
2642 // Sanitize dataset name
2643 data->SetName(sanitizeName(data->GetName()).c_str());
2644 for (auto *obj : *data->get()) {
2645 obj->SetName(sanitizeName(obj->GetName()).c_str());
2646 }
2647 }
2648 for (auto *snshObj : tmpWS.getSnapshots()) {
2649 // Snapshots are stored as TObject*, but really they are RooArgSet*
2650 auto *snsh = dynamic_cast<RooArgSet *>(snshObj);
2651 if (!snsh) {
2652 std::cerr << "Warning: found snapshot that is not a RooArgSet, skipping\n";
2653 continue;
2654 }
2655
2656 // Sanitize snapshot name
2657 if (!isValidName(snsh->GetName())) {
2658 snsh->setName(sanitizeName(snsh->GetName()).c_str());
2659 }
2660
2661 // Sanitize the variables inside the snapshot
2662 for (auto *arg : *snsh) {
2663 if (!isValidName(arg->GetName())) {
2664 arg->SetName(sanitizeName(arg->GetName()).c_str());
2665 }
2666 }
2667 }
2668
2669 // Generic objects (ModelConfigs, attributes, etc.)
2670 for (auto *obj : tmpWS.allGenericObjects()) {
2671 if (!isValidName(obj->GetName())) {
2672 if (auto *named = dynamic_cast<TNamed *>(obj)) {
2673 named->SetName(sanitizeName(named->GetName()).c_str());
2674 } else {
2675 std::cerr << "Warning: object " << obj->GetName() << " is not TNamed, cannot rename.\n";
2676 }
2677 }
2678
2679 if (auto *mc = dynamic_cast<RooFit::ModelConfig *>(obj)) {
2680 // Sanitize ModelConfig name
2681 if (!isValidName(mc->GetName())) {
2682 mc->SetName(sanitizeName(mc->GetName()).c_str());
2683 }
2684
2685 // Sanitize the sets inside ModelConfig
2686 for (auto *obs : mc->GetObservables()->get()) {
2687 if (obs) {
2688 obs->SetName(sanitizeName(obs->GetName()).c_str());
2689 }
2690 }
2691 for (auto *poi : mc->GetParametersOfInterest()->get()) {
2692 if (poi) {
2693 poi->SetName(sanitizeName(poi->GetName()).c_str());
2694 }
2695 }
2696 for (auto *nuis : mc->GetNuisanceParameters()->get()) {
2697 if (nuis) {
2698 nuis->SetName(sanitizeName(nuis->GetName()).c_str());
2699 }
2700 }
2701 for (auto *glob : mc->GetGlobalObservables()->get()) {
2702 if (glob) {
2703 glob->SetName(sanitizeName(glob->GetName()).c_str());
2704 }
2705 }
2706 }
2707 }
2708 std::string wsName = std::string{ws.GetName()} + "_sanitized";
2709 RooWorkspace newWS = cleanWS(tmpWS, false);
2710 newWS.SetName(wsName.c_str());
2711
2712 return newWS;
2713}
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
#define e(i)
Definition RSha256.hxx:103
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 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:148
#define gROOT
Definition TROOT.h:417
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...
RooFit::OwningPtr< RooArgSet > getObservables(const RooArgSet &set, bool valueOnly=true) const
Given a set of possible observables, return the observables that this PDF depends on.
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:56
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
void fill_seq(Collection const &coll)
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
static std::unique_ptr< JSONTree > create()
void writeJSON(RooFit::Detail::JSONNode &) const
Definition Domains.cxx:144
When using RooFit, statistical models can be conveniently handled and stored as a RooWorkspace.
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 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)
bool importJSON(std::string const &filename)
Imports a JSON file from the given filename to the workspace.
void exportVariables(const RooArgSet &allElems, RooFit::Detail::JSONNode &n, bool storeConstant, bool storeBins)
Export variables from the workspace to a JSONNode.
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
std::unordered_map< std::string, RooFit::Detail::JSONNode const * > _functionsByName
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.
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)
void importVariable(const RooFit::Detail::JSONNode &p)
Import a variable from the JSONNode into the workspace.
bool exportYML(std::string const &fileName)
Export the workspace to YML format and write to the specified file.
void exportVariable(const RooAbsArg *v, RooFit::Detail::JSONNode &n, bool storeConstant, bool storeBins)
Export a variable from the workspace to a JSONNode.
void importFunction(const RooFit::Detail::JSONNode &p, bool importAllDependants)
Import a function from the JSONNode into the workspace.
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)
std::unordered_map< std::string, RooFit::Detail::JSONNode const * > _distributionsByName
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.
RooArgSet allResolutionModels() const
Return set with all resolution model objects.
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:42
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:2385
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:60
ExportMap & exporters()
Definition JSONIO.cxx:82
ImportExpressionMap & importExpressions()
Definition JSONIO.cxx:109
ExportKeysMap & exportKeys()
Definition JSONIO.cxx:116
TLine l
Definition textangle.C:4