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