Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
JSONFactories_RooFitCore.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
14
15#include <RooAbsCachedPdf.h>
16#include <RooAddPdf.h>
17#include <RooAddModel.h>
18#include <RooBinning.h>
19#include <RooBinSamplingPdf.h>
20#include <RooBinWidthFunction.h>
21#include <RooCategory.h>
22#include <RooDataHist.h>
23#include <RooDecay.h>
24#include <RooDerivative.h>
25#include <RooExponential.h>
26#include <RooExtendPdf.h>
27#include <RooFFTConvPdf.h>
29#include <RooFitHS3/JSONIO.h>
30#include <RooFormulaVar.h>
31#include <RooGenericPdf.h>
32#include <RooHistFunc.h>
33#include <RooHistPdf.h>
34#include <RooLegacyExpPoly.h>
35#include <RooLognormal.h>
36#include <RooMultiVarGaussian.h>
38#include <RooPoisson.h>
39#include <RooPolynomial.h>
40#include <RooPolyVar.h>
41#include <RooRealSumFunc.h>
42#include <RooRealSumPdf.h>
43#include <RooRealVar.h>
44#include <RooResolutionModel.h>
45#include <RooTFnBinding.h>
46#include <RooTruthModel.h>
47#include <RooGaussModel.h>
48#include <RooWorkspace.h>
49#include <RooRealIntegral.h>
50#include <RooSpline.h>
51#include <TSpline.h>
52
53#include <TF1.h>
54#include <TH1.h>
55
56#include "JSONIOUtils.h"
57
58#include "static_execute.h"
59
60#include <algorithm>
61#include <cctype>
62
64
65///////////////////////////////////////////////////////////////////////////////////////////////////////
66// individually implemented importers
67///////////////////////////////////////////////////////////////////////////////////////////////////////
68
69namespace {
70/**
71 * Extracts arguments from a mathematical expression.
72 *
73 * This function takes a string representing a mathematical
74 * expression and extracts the arguments from it. The arguments are
75 * defined as sequences of characters that do not contain digits,
76 * spaces, or parentheses, and that start with a letter. Function
77 * calls such as "exp( ... )", identified as being followed by an
78 * opening parenthesis, are not treated as arguments. The extracted
79 * arguments are returned as a vector of strings.
80 *
81 * @param expr A string representing a mathematical expression.
82 * @return A set of unique strings representing the extracted arguments.
83 */
84std::set<std::string> extractArguments(std::string expr)
85{
86 // Get rid of whitespaces
87 expr.erase(std::remove_if(expr.begin(), expr.end(), [](unsigned char c) { return std::isspace(c); }), expr.end());
88
89 std::set<std::string> arguments;
90 size_t startidx = expr.size();
91 for (size_t i = 0; i < expr.size(); ++i) {
92 if (startidx >= expr.size()) {
93 if (isalpha(expr[i])) {
94 startidx = i;
95 // check this character is not part of scientific notation, e.g. 2e-5
97 // if it is, we ignore this character
98 startidx = expr.size();
99 }
100 }
101 } else {
102 if (!isdigit(expr[i]) && !isalpha(expr[i]) && expr[i] != '_') {
103 if (expr[i] == '(') {
104 startidx = expr.size();
105 continue;
106 }
107 std::string arg(expr.substr(startidx, i - startidx));
108 startidx = expr.size();
109 arguments.insert(arg);
110 }
111 }
112 }
113 if (startidx < expr.size()) {
114 arguments.insert(expr.substr(startidx));
115 }
116 return arguments;
117}
118
119template <class RooArg_t>
120class RooFormulaArgFactory : public RooFit::JSONIO::Importer {
121public:
122 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
123 {
124 std::string name(RooJSONFactoryWSTool::name(p));
125 if (!p.has_child("expression")) {
126 RooJSONFactoryWSTool::error("no expression given for '" + name + "'");
127 }
128 TString formula(p["expression"].val());
129 RooArgList dependents;
130 for (const auto &d : extractArguments(formula.Data())) {
131 dependents.add(*tool->request<RooAbsReal>(d, name));
132 }
133 tool->wsImport(RooArg_t{name.c_str(), formula, dependents});
134 return true;
135 }
136};
137
138class RooAddPdfFactory : public RooFit::JSONIO::Importer {
139public:
140 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
141 {
142 std::string name(RooJSONFactoryWSTool::name(p));
143 if (!tool->requestArgList<RooAbsReal>(p, "coefficients").empty()) {
144 tool->wsEmplace<RooAddPdf>(name, tool->requestArgList<RooAbsPdf>(p, "summands"),
145 tool->requestArgList<RooAbsReal>(p, "coefficients"));
146 return true;
147 }
148 tool->wsEmplace<RooAddPdf>(name, tool->requestArgList<RooAbsPdf>(p, "summands"));
149 return true;
150 }
151};
152
153class RooAddModelFactory : public RooFit::JSONIO::Importer {
154public:
155 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
156 {
157 std::string name(RooJSONFactoryWSTool::name(p));
158 tool->wsEmplace<RooAddModel>(name, tool->requestArgList<RooAbsPdf>(p, "summands"),
159 tool->requestArgList<RooAbsReal>(p, "coefficients"));
160 return true;
161 }
162};
163
164class RooBinWidthFunctionFactory : public RooFit::JSONIO::Importer {
165public:
166 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
167 {
168 std::string name(RooJSONFactoryWSTool::name(p));
169 RooHistFunc *hf = static_cast<RooHistFunc *>(tool->request<RooAbsReal>(p["histogram"].val(), name));
170 tool->wsEmplace<RooBinWidthFunction>(name, *hf, p["divideByBinWidth"].val_bool());
171 return true;
172 }
173};
174
175class RooBinSamplingPdfFactory : public RooFit::JSONIO::Importer {
176public:
177 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
178 {
179 std::string name(RooJSONFactoryWSTool::name(p));
180
181 RooAbsPdf *pdf = tool->requestArg<RooAbsPdf>(p, "pdf");
182 RooRealVar *obs = tool->requestArg<RooRealVar>(p, "observable");
183
184 if (!pdf->dependsOn(*obs)) {
185 RooJSONFactoryWSTool::error(std::string("pdf '") + pdf->GetName() + "' does not depend on observable '" +
186 obs->GetName() + "' as indicated by parent RooBinSamplingPdf '" + name +
187 "', please check!");
188 }
189
190 if (!p.has_child("epsilon")) {
191 RooJSONFactoryWSTool::error("no epsilon given in '" + name + "'");
192 }
193 double epsilon(p["epsilon"].val_double());
194
195 tool->wsEmplace<RooBinSamplingPdf>(name, *obs, *pdf, epsilon);
196
197 return true;
198 }
199};
200
201class RooRealSumPdfFactory : public RooFit::JSONIO::Importer {
202public:
203 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
204 {
205 std::string name(RooJSONFactoryWSTool::name(p));
206
207 bool extended = false;
208 if (p.has_child("extended") && p["extended"].val_bool()) {
209 extended = true;
210 }
211 tool->wsEmplace<RooRealSumPdf>(name, tool->requestArgList<RooAbsReal>(p, "samples"),
212 tool->requestArgList<RooAbsReal>(p, "coefficients"), extended);
213 return true;
214 }
215};
216
217class RooRealSumFuncFactory : public RooFit::JSONIO::Importer {
218public:
219 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
220 {
221 std::string name(RooJSONFactoryWSTool::name(p));
222 tool->wsEmplace<RooRealSumFunc>(name, tool->requestArgList<RooAbsReal>(p, "samples"),
223 tool->requestArgList<RooAbsReal>(p, "coefficients"));
224 return true;
225 }
226};
227template <class RooArg_t>
228class RooPolynomialFactory : public RooFit::JSONIO::Importer {
229public:
230 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
231 {
232 std::string name(RooJSONFactoryWSTool::name(p));
233 if (!p.has_child("coefficients")) {
234 RooJSONFactoryWSTool::error("no coefficients given in '" + name + "'");
235 }
236 RooAbsReal *x = tool->requestArg<RooAbsReal>(p, "x");
237 RooArgList coefs;
238 int order = 0;
239 int lowestOrder = 0;
240 for (const auto &coef : p["coefficients"].children()) {
241 // As long as the coefficients match the default coefficients in
242 // RooFit, we don't have to instantiate RooFit objects but can
243 // increase the lowestOrder flag.
244 if (order == 0 && coef.val() == "1.0") {
245 ++lowestOrder;
246 } else if (coefs.empty() && coef.val() == "0.0") {
247 ++lowestOrder;
248 } else {
249 coefs.add(*tool->request<RooAbsReal>(coef.val(), name));
250 }
251 ++order;
252 }
253
254 tool->wsEmplace<RooArg_t>(name, *x, coefs, lowestOrder);
255 return true;
256 }
257};
258
259class RooPoissonFactory : public RooFit::JSONIO::Importer {
260public:
261 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
262 {
263 std::string name(RooJSONFactoryWSTool::name(p));
264 RooAbsReal *x = tool->requestArg<RooAbsReal>(p, "x");
265 RooAbsReal *mean = tool->requestArg<RooAbsReal>(p, "mean");
266 tool->wsEmplace<RooPoisson>(name, *x, *mean, !p["integer"].val_bool());
267 return true;
268 }
269};
270
271class RooDecayFactory : public RooFit::JSONIO::Importer {
272public:
273 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
274 {
275 std::string name(RooJSONFactoryWSTool::name(p));
276 RooRealVar *t = tool->requestArg<RooRealVar>(p, "t");
277 RooAbsReal *tau = tool->requestArg<RooAbsReal>(p, "tau");
278 RooResolutionModel *model = dynamic_cast<RooResolutionModel *>(tool->requestArg<RooAbsPdf>(p, "resolutionModel"));
279 RooDecay::DecayType decayType = static_cast<RooDecay::DecayType>(p["decayType"].val_int());
280 tool->wsEmplace<RooDecay>(name, *t, *tau, *model, decayType);
281 return true;
282 }
283};
284
285class RooTruthModelFactory : public RooFit::JSONIO::Importer {
286public:
287 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
288 {
289 std::string name(RooJSONFactoryWSTool::name(p));
290 RooRealVar *x = tool->requestArg<RooRealVar>(p, "x");
291 tool->wsEmplace<RooTruthModel>(name, *x);
292 return true;
293 }
294};
295
296class RooGaussModelFactory : public RooFit::JSONIO::Importer {
297public:
298 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
299 {
300 std::string name(RooJSONFactoryWSTool::name(p));
301 RooRealVar *x = tool->requestArg<RooRealVar>(p, "x");
302 RooRealVar *mean = tool->requestArg<RooRealVar>(p, "mean");
303 RooRealVar *sigma = tool->requestArg<RooRealVar>(p, "sigma");
304 tool->wsEmplace<RooGaussModel>(name, *x, *mean, *sigma);
305 return true;
306 }
307};
308
309class RooRealIntegralFactory : public RooFit::JSONIO::Importer {
310public:
311 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
312 {
313 std::string name(RooJSONFactoryWSTool::name(p));
314 RooAbsReal *func = tool->requestArg<RooAbsReal>(p, "integrand");
315 auto vars = tool->requestArgList<RooAbsReal>(p, "variables");
317 RooArgSet const *normSetPtr = nullptr;
318 if (p.has_child("normalization")) {
319 normSet.add(tool->requestArgSet<RooAbsReal>(p, "normalization"));
321 }
322 std::string domain;
323 bool hasDomain = p.has_child("domain");
324 if (hasDomain) {
325 domain = p["domain"].val();
326 }
327 // todo: at some point, take care of integrator configurations
328 tool->wsEmplace<RooRealIntegral>(name, *func, vars, normSetPtr, static_cast<RooNumIntConfig *>(nullptr),
329 hasDomain ? domain.c_str() : nullptr);
330 return true;
331 }
332};
333
334class RooDerivativeFactory : public RooFit::JSONIO::Importer {
335public:
336 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
337 {
338 std::string name(RooJSONFactoryWSTool::name(p));
339 RooAbsReal *func = tool->requestArg<RooAbsReal>(p, "function");
340 RooRealVar *x = tool->requestArg<RooRealVar>(p, "x");
341 Int_t order = p["order"].val_int();
342 double eps = p["eps"].val_double();
343 if (p.has_child("normalization")) {
345 normSet.add(tool->requestArgSet<RooAbsReal>(p, "normalization"));
346 tool->wsEmplace<RooDerivative>(name, *func, *x, normSet, order, eps);
347 return true;
348 }
349 tool->wsEmplace<RooDerivative>(name, *func, *x, order, eps);
350 return true;
351 }
352};
353
354class RooFFTConvPdfFactory : public RooFit::JSONIO::Importer {
355public:
356 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
357 {
358 std::string name(RooJSONFactoryWSTool::name(p));
359 RooRealVar *convVar = tool->requestArg<RooRealVar>(p, "conv_var");
360 Int_t order = p["ipOrder"].val_int();
361 RooAbsPdf *pdf1 = tool->requestArg<RooAbsPdf>(p, "pdf1");
362 RooAbsPdf *pdf2 = tool->requestArg<RooAbsPdf>(p, "pdf2");
363 if (p.has_child("conv_func")) {
364 RooAbsReal *convFunc = tool->requestArg<RooAbsReal>(p, "conv_func");
365 tool->wsEmplace<RooFFTConvPdf>(name, *convFunc, *convVar, *pdf1, *pdf2, order);
366 return true;
367 }
368 tool->wsEmplace<RooFFTConvPdf>(name, *convVar, *pdf1, *pdf2, order);
369 return true;
370 }
371};
372
373class RooExtendPdfFactory : public RooFit::JSONIO::Importer {
374public:
375 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
376 {
377 std::string name(RooJSONFactoryWSTool::name(p));
378 RooAbsPdf *pdf = tool->requestArg<RooAbsPdf>(p, "pdf");
379 RooAbsReal *norm = tool->requestArg<RooAbsReal>(p, "norm");
380 if (p.has_child("range")) {
381 std::string rangeName = p["range"].val();
382 tool->wsEmplace<RooExtendPdf>(name, *pdf, *norm, rangeName.c_str());
383 return true;
384 }
385 tool->wsEmplace<RooExtendPdf>(name, *pdf, *norm);
386 return true;
387 }
388};
389
390class RooLogNormalFactory : public RooFit::JSONIO::Importer {
391public:
392 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
393 {
394 std::string name(RooJSONFactoryWSTool::name(p));
395 RooAbsReal *x = tool->requestArg<RooAbsReal>(p, "x");
396
397 // Same mechanism to undo the parameter transformation as in the
398 // RooExponentialFactory (see comments in that class for more info).
399 const std::string muName = p["mu"].val();
400 const std::string sigmaName = p["sigma"].val();
401 const bool isTransformed = endsWith(muName, "_lognormal_log");
402 const std::string suffixToRemove = isTransformed ? "_lognormal_log" : "";
405
406 tool->wsEmplace<RooLognormal>(name, *x, *mu, *sigma, !isTransformed);
407
408 return true;
409 }
410};
411
412class RooExponentialFactory : public RooFit::JSONIO::Importer {
413public:
414 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
415 {
416 std::string name(RooJSONFactoryWSTool::name(p));
417 RooAbsReal *x = tool->requestArg<RooAbsReal>(p, "x");
418
419 // If the parameter name ends with the "_exponential_inverted" suffix,
420 // this means that it was exported from a RooFit object where the
421 // parameter first needed to be transformed on export to match the HS3
422 // specification. But when re-importing such a parameter, we can simply
423 // skip the transformation and use the original RooFit parameter without
424 // the suffix.
425 //
426 // A concrete example: take the following RooFit pdf in the factory language:
427 //
428 // "Exponential::exponential_1(x[0, 10], c[-0.1])"
429 //
430 // It defines en exponential exp(c * x). However, in HS3 the exponential
431 // is defined as exp(-c * x), to RooFit would export these dictionaries
432 // to the JSON:
433 //
434 // {
435 // "name": "exponential_1", // HS3 exponential_dist with transformed parameter
436 // "type": "exponential_dist",
437 // "x": "x",
438 // "c": "c_exponential_inverted"
439 // },
440 // {
441 // "name": "c_exponential_inverted", // transformation function created on-the-fly on export
442 // "type": "generic_function",
443 // "expression": "-c"
444 // }
445 //
446 // On import, we can directly take the non-transformed parameter, which is
447 // we check for the suffix and optionally remove it from the requested
448 // name next:
449
450 const std::string constParamName = p["c"].val();
451 const bool isInverted = endsWith(constParamName, "_exponential_inverted");
452 const std::string suffixToRemove = isInverted ? "_exponential_inverted" : "";
454
455 tool->wsEmplace<RooExponential>(name, *x, *c, !isInverted);
456
457 return true;
458 }
459};
460
461class RooLegacyExpPolyFactory : public RooFit::JSONIO::Importer {
462public:
463 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
464 {
465 std::string name(RooJSONFactoryWSTool::name(p));
466 if (!p.has_child("coefficients")) {
467 RooJSONFactoryWSTool::error("no coefficients given in '" + name + "'");
468 }
469 RooAbsReal *x = tool->requestArg<RooAbsReal>(p, "x");
470 RooArgList coefs;
471 int order = 0;
472 int lowestOrder = 0;
473 for (const auto &coef : p["coefficients"].children()) {
474 // As long as the coefficients match the default coefficients in
475 // RooFit, we don't have to instantiate RooFit objects but can
476 // increase the lowestOrder flag.
477 if (order == 0 && coef.val() == "1.0") {
478 ++lowestOrder;
479 } else if (coefs.empty() && coef.val() == "0.0") {
480 ++lowestOrder;
481 } else {
482 coefs.add(*tool->request<RooAbsReal>(coef.val(), name));
483 }
484 ++order;
485 }
486
487 tool->wsEmplace<RooLegacyExpPoly>(name, *x, coefs, lowestOrder);
488 return true;
489 }
490};
491
492class RooMultiVarGaussianFactory : public RooFit::JSONIO::Importer {
493public:
494 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
495 {
496 std::string name(RooJSONFactoryWSTool::name(p));
497 bool has_cov = p.has_child("covariances");
498 bool has_corr = p.has_child("correlations") && p.has_child("standard_deviations");
499 if (!has_cov && !has_corr) {
500 RooJSONFactoryWSTool::error("no covariances or correlations+standard_deviations given in '" + name + "'");
501 }
502
504
505 if (has_cov) {
506 int n = p["covariances"].num_children();
507 int i = 0;
508 covmat.ResizeTo(n, n);
509 for (const auto &row : p["covariances"].children()) {
510 int j = 0;
511 for (const auto &val : row.children()) {
512 covmat(i, j) = val.val_double();
513 ++j;
514 }
515 ++i;
516 }
517 } else {
518 std::vector<double> variances;
519 for (const auto &v : p["standard_deviations"].children()) {
520 variances.push_back(v.val_double());
521 }
522 covmat.ResizeTo(variances.size(), variances.size());
523 int i = 0;
524 for (const auto &row : p["correlations"].children()) {
525 int j = 0;
526 for (const auto &val : row.children()) {
527 covmat(i, j) = val.val_double() * variances[i] * variances[j];
528 ++j;
529 }
530 ++i;
531 }
532 }
533 tool->wsEmplace<RooMultiVarGaussian>(name, tool->requestArgList<RooAbsReal>(p, "x"),
534 tool->requestArgList<RooAbsReal>(p, "mean"), covmat);
535 return true;
536 }
537};
538
539class ParamHistFuncFactory : public RooFit::JSONIO::Importer {
540public:
541 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
542 {
543 std::string name(RooJSONFactoryWSTool::name(p));
544 RooArgList varList = tool->requestArgList<RooRealVar>(p, "variables");
545 if (!p.has_child("axes")) {
546 std::stringstream ss;
547 ss << "No axes given in '" << name << "'"
548 << ". Using default binning (uniform; nbins=100). If needed, export the Workspace to JSON with a newer "
549 << "Root version that supports custom ParamHistFunc binnings(>=6.38.00)." << std::endl;
551 tool->wsEmplace<ParamHistFunc>(name, varList, tool->requestArgList<RooAbsReal>(p, "parameters"));
552 return true;
553 }
554 tool->wsEmplace<ParamHistFunc>(name, readBinning(p, varList), tool->requestArgList<RooAbsReal>(p, "parameters"));
555 return true;
556 }
557
558private:
560 {
561 // Temporary map from variable name → RooRealVar
562 std::map<std::string, std::unique_ptr<RooRealVar>> varMap;
563
564 // Build variables from JSON
565 for (const JSONNode &node : topNode["axes"].children()) {
566 const std::string name = node["name"].val();
567 std::unique_ptr<RooRealVar> obs;
568
569 if (node.has_child("edges")) {
570 std::vector<double> edges;
571 for (const auto &bound : node["edges"].children()) {
572 edges.push_back(bound.val_double());
573 }
574 obs = std::make_unique<RooRealVar>(name.c_str(), name.c_str(), edges.front(), edges.back());
575 RooBinning bins(obs->getMin(), obs->getMax());
576 for (auto b : edges)
577 bins.addBoundary(b);
578 obs->setBinning(bins);
579 } else {
580 obs = std::make_unique<RooRealVar>(name.c_str(), name.c_str(), node["min"].val_double(),
581 node["max"].val_double());
582 obs->setBins(node["nbins"].val_int());
583 }
584
585 varMap[name] = std::move(obs);
586 }
587
588 // Now build the final list following the order in varList
589 RooArgList vars;
590 for (int i = 0; i < varList.getSize(); ++i) {
591 const auto *refVar = dynamic_cast<RooRealVar *>(varList.at(i));
592 if (!refVar)
593 continue;
594
595 auto it = varMap.find(refVar->GetName());
596 if (it != varMap.end()) {
597 vars.addOwned(std::move(it->second)); // preserve ownership
598 }
599 }
600 return vars;
601 }
602};
603
604class RooSplineFactory : public RooFit::JSONIO::Importer {
605public:
606 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
607 {
608 const std::string name(RooJSONFactoryWSTool::name(p));
609
610 // Mandatory fields
611 if (!p.has_child("x")) {
612 RooJSONFactoryWSTool::error("no x given in '" + name + "'");
613 }
614 if (!p.has_child("x0") || !p.has_child("y0")) {
615 RooJSONFactoryWSTool::error("no x0/y0 given in '" + name + "'");
616 }
617
618 RooAbsReal *x = tool->requestArg<RooAbsReal>(p, "x");
619
620 // Optional fields (defaults follow RooSpline ctor defaults)
621 std::string algo = p.has_child("interpolation") ? p["interpolation"].val() : "poly3";
622 int order = 0;
623 if (algo == "poly3")
624 order = 3;
625 else if (algo == "poly5")
626 order = 5;
627 else {
628 RooJSONFactoryWSTool::error("unsupported algo '" + algo + "' for RooSpline in '" + name +
629 "': allowed are 'poly3' and 'poly5'");
630 }
631 const bool logx = p.has_child("logx") ? p["logx"].val_bool() : false;
632 const bool logy = p.has_child("logy") ? p["logy"].val_bool() : false;
633
634 // Read knots
635 std::vector<double> x0;
636 std::vector<double> y0;
637 x0.reserve(p["x0"].num_children());
638 y0.reserve(p["y0"].num_children());
639
640 for (const auto &v : p["x0"].children())
641 x0.push_back(v.val_double());
642 for (const auto &v : p["y0"].children())
643 y0.push_back(v.val_double());
644
645 if (x0.size() != y0.size()) {
646 RooJSONFactoryWSTool::error("x0/y0 size mismatch in '" + name + "': x0 has " + std::to_string(x0.size()) +
647 ", y0 has " + std::to_string(y0.size()));
648 }
649 if (x0.size() < 2) {
650 RooJSONFactoryWSTool::error("need at least 2 knots in '" + name + "'");
651 }
652
653 // Construct RooSpline(name,title, x, x0, y0, order, logx, logy)
654 tool->wsEmplace<::RooSpline>(name.c_str(), *x, std::span<const double>(x0.data(), x0.size()),
655 std::span<const double>(y0.data(), y0.size()), order, logx, logy);
656
657 return true;
658 }
659};
660
661///////////////////////////////////////////////////////////////////////////////////////////////////////
662// specialized exporter implementations
663///////////////////////////////////////////////////////////////////////////////////////////////////////
664template <class RooArg_t>
665class RooAddPdfStreamer : public RooFit::JSONIO::Exporter {
666public:
667 std::string const &key() const override;
668 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
669 {
670 const RooArg_t *pdf = static_cast<const RooArg_t *>(func);
671 elem["type"] << key();
672 RooJSONFactoryWSTool::fillSeq(elem["summands"], pdf->pdfList());
673 RooJSONFactoryWSTool::fillSeq(elem["coefficients"], pdf->coefList());
674 elem["extended"] << (pdf->extendMode() != RooArg_t::CanNotBeExtended);
675 return true;
676 }
677};
678
679class RooRealSumPdfStreamer : public RooFit::JSONIO::Exporter {
680public:
681 std::string const &key() const override;
682 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
683 {
684 const RooRealSumPdf *pdf = static_cast<const RooRealSumPdf *>(func);
685 elem["type"] << key();
686 RooJSONFactoryWSTool::fillSeq(elem["samples"], pdf->funcList());
687 RooJSONFactoryWSTool::fillSeq(elem["coefficients"], pdf->coefList());
688 elem["extended"] << (pdf->extendMode() != RooAbsPdf::CanNotBeExtended);
689 return true;
690 }
691};
692
693class RooRealSumFuncStreamer : public RooFit::JSONIO::Exporter {
694public:
695 std::string const &key() const override;
696 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
697 {
698 const RooRealSumFunc *pdf = static_cast<const RooRealSumFunc *>(func);
699 elem["type"] << key();
700 RooJSONFactoryWSTool::fillSeq(elem["samples"], pdf->funcList());
701 RooJSONFactoryWSTool::fillSeq(elem["coefficients"], pdf->coefList());
702 return true;
703 }
704};
705
706class RooHistFuncStreamer : public RooFit::JSONIO::Exporter {
707public:
708 std::string const &key() const override;
709 bool exportObject(RooJSONFactoryWSTool *tool, const RooAbsArg *func, JSONNode &elem) const override
710 {
711 const RooHistFunc *hf = static_cast<const RooHistFunc *>(func);
712 elem["type"] << key();
713 RooDataHist const &dh = hf->dataHist();
714 tool->exportHisto(*dh.get(), dh.numEntries(), dh.weightArray(), elem["data"].set_map());
715 return true;
716 }
717};
718
719class RooHistFuncFactory : public RooFit::JSONIO::Importer {
720public:
721 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
722 {
723 std::string name(RooJSONFactoryWSTool::name(p));
724 if (!p.has_child("data")) {
725 RooJSONFactoryWSTool::error("function '" + name + "' is of histogram type, but does not define a 'data' key");
726 }
727 std::unique_ptr<RooDataHist> dataHist =
729 tool->wsEmplace<RooHistFunc>(name, *dataHist->get(), *dataHist);
730 return true;
731 }
732};
733
734class RooHistPdfStreamer : public RooFit::JSONIO::Exporter {
735public:
736 std::string const &key() const override;
737 bool exportObject(RooJSONFactoryWSTool *tool, const RooAbsArg *func, JSONNode &elem) const override
738 {
739 const RooHistPdf *hf = static_cast<const RooHistPdf *>(func);
740 elem["type"] << key();
741 RooDataHist const &dh = hf->dataHist();
742 tool->exportHisto(*dh.get(), dh.numEntries(), dh.weightArray(), elem["data"].set_map());
743 return true;
744 }
745};
746
747class RooHistPdfFactory : public RooFit::JSONIO::Importer {
748public:
749 bool importArg(RooJSONFactoryWSTool *tool, const JSONNode &p) const override
750 {
751 std::string name(RooJSONFactoryWSTool::name(p));
752 if (!p.has_child("data")) {
753 RooJSONFactoryWSTool::error("function '" + name + "' is of histogram type, but does not define a 'data' key");
754 }
755 std::unique_ptr<RooDataHist> dataHist =
757 tool->wsEmplace<RooHistPdf>(name, *dataHist->get(), *dataHist);
758 return true;
759 }
760};
761
762class RooBinSamplingPdfStreamer : public RooFit::JSONIO::Exporter {
763public:
764 std::string const &key() const override;
765 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
766 {
767 const RooBinSamplingPdf *pdf = static_cast<const RooBinSamplingPdf *>(func);
768 elem["type"] << key();
769 elem["pdf"] << pdf->pdf().GetName();
770 elem["observable"] << pdf->observable().GetName();
771 elem["epsilon"] << pdf->epsilon();
772 return true;
773 }
774};
775
776class RooBinWidthFunctionStreamer : public RooFit::JSONIO::Exporter {
777public:
778 std::string const &key() const override;
779 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
780 {
781 const RooBinWidthFunction *pdf = static_cast<const RooBinWidthFunction *>(func);
782 elem["type"] << key();
783 elem["histogram"] << pdf->histFunc().GetName();
784 elem["divideByBinWidth"] << pdf->divideByBinWidth();
785 return true;
786 }
787};
788
789template <class RooArg_t>
790class RooFormulaArgStreamer : public RooFit::JSONIO::Exporter {
791public:
792 std::string const &key() const override;
793 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
794 {
795 const RooArg_t *pdf = static_cast<const RooArg_t *>(func);
796 elem["type"] << key();
797 TString expression(pdf->expression());
798 cleanExpression(expression);
799 // If the tokens follow the "x[#]" convention, the square braces enclosing each number
800 // ensures that there is a unique mapping between the token and parameter name
801 // If the tokens follow the "@#" convention, the numbers are not enclosed by braces.
802 // So there may be tokens with numbers whose lower place value forms a subset string of ones with a higher place
803 // value, e.g. "@1" is a subset of "@10". So the names of these parameters must be applied descending from the
804 // highest place value in order to ensure each parameter name is uniquely applied to its token.
805 for (size_t idx = pdf->nParameters(); idx--;) {
806 const RooAbsArg *par = pdf->getParameter(idx);
807 expression.ReplaceAll(("x[" + std::to_string(idx) + "]").c_str(), par->GetName());
808 expression.ReplaceAll(("@" + std::to_string(idx)).c_str(), par->GetName());
809 }
810 elem["expression"] << expression.Data();
811 return true;
812 }
813
814private:
815 void cleanExpression(TString &expr) const
816 {
817 expr.ReplaceAll("TMath::Exp", "exp");
818 expr.ReplaceAll("TMath::Min", "min");
819 expr.ReplaceAll("TMath::Max", "max");
820 expr.ReplaceAll("TMath::Log", "log");
821 expr.ReplaceAll("TMath::Cos", "cos");
822 expr.ReplaceAll("TMath::Sin", "sin");
823 expr.ReplaceAll("TMath::Sqrt", "sqrt");
824 expr.ReplaceAll("TMath::Power", "pow");
825 expr.ReplaceAll("TMath::Erf", "erf");
826 }
827};
828template <class RooArg_t>
829class RooPolynomialStreamer : public RooFit::JSONIO::Exporter {
830public:
831 std::string const &key() const override;
832 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
833 {
834 auto *pdf = static_cast<const RooArg_t *>(func);
835 elem["type"] << key();
836 elem["x"] << pdf->x().GetName();
837 auto &coefs = elem["coefficients"].set_seq();
838 // Write out the default coefficient that RooFit uses for the lower
839 // orders before the order of the first coefficient. Like this, the
840 // output is more self-documenting.
841 for (int i = 0; i < pdf->lowestOrder(); ++i) {
842 coefs.append_child() << (i == 0 ? "1.0" : "0.0");
843 }
844 for (const auto &coef : pdf->coefList()) {
845 coefs.append_child() << coef->GetName();
846 }
847 return true;
848 }
849};
850
851class RooLegacyExpPolyStreamer : public RooFit::JSONIO::Exporter {
852public:
853 std::string const &key() const override;
854 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
855 {
856 auto *pdf = static_cast<const RooLegacyExpPoly *>(func);
857 elem["type"] << key();
858 elem["x"] << pdf->x().GetName();
859 auto &coefs = elem["coefficients"].set_seq();
860 // Write out the default coefficient that RooFit uses for the lower
861 // orders before the order of the first coefficient. Like this, the
862 // output is more self-documenting.
863 for (int i = 0; i < pdf->lowestOrder(); ++i) {
864 coefs.append_child() << (i == 0 ? "1.0" : "0.0");
865 }
866 for (const auto &coef : pdf->coefList()) {
867 coefs.append_child() << coef->GetName();
868 }
869 return true;
870 }
871};
872
873class RooPoissonStreamer : public RooFit::JSONIO::Exporter {
874public:
875 std::string const &key() const override;
876 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
877 {
878 auto *pdf = static_cast<const RooPoisson *>(func);
879 elem["type"] << key();
880 elem["x"] << pdf->getX().GetName();
881 elem["mean"] << pdf->getMean().GetName();
882 elem["integer"] << !pdf->getNoRounding();
883 return true;
884 }
885};
886
887class RooDecayStreamer : public RooFit::JSONIO::Exporter {
888public:
889 std::string const &key() const override;
890 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
891 {
892 auto *pdf = static_cast<const RooDecay *>(func);
893 elem["type"] << key();
894 elem["t"] << pdf->getT().GetName();
895 elem["tau"] << pdf->getTau().GetName();
896 elem["resolutionModel"] << pdf->getModel().GetName();
897 elem["decayType"] << pdf->getDecayType();
898
899 return true;
900 }
901};
902
903class RooTruthModelStreamer : public RooFit::JSONIO::Exporter {
904public:
905 std::string const &key() const override;
906 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
907 {
908 auto *pdf = static_cast<const RooTruthModel *>(func);
909 elem["type"] << key();
910 elem["x"] << pdf->convVar().GetName();
911
912 return true;
913 }
914};
915
916class RooGaussModelStreamer : public RooFit::JSONIO::Exporter {
917public:
918 std::string const &key() const override;
919 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
920 {
921 auto *pdf = static_cast<const RooGaussModel *>(func);
922 elem["type"] << key();
923 elem["x"] << pdf->convVar().GetName();
924 elem["mean"] << pdf->getMean().GetName();
925 elem["sigma"] << pdf->getSigma().GetName();
926 return true;
927 }
928};
929
930class RooLogNormalStreamer : public RooFit::JSONIO::Exporter {
931public:
932 std::string const &key() const override;
933 bool exportObject(RooJSONFactoryWSTool *tool, const RooAbsArg *func, JSONNode &elem) const override
934 {
935 auto *pdf = static_cast<const RooLognormal *>(func);
936
937 elem["type"] << key();
938 elem["x"] << pdf->getX().GetName();
939
940 auto &m0 = pdf->getMedian();
941 auto &k = pdf->getShapeK();
942
943 if (pdf->useStandardParametrization()) {
944 elem["mu"] << m0.GetName();
945 elem["sigma"] << k.GetName();
946 } else {
947 elem["mu"] << tool->exportTransformed(&m0, "_lognormal_log", "log(%s)");
948 elem["sigma"] << tool->exportTransformed(&k, "_lognormal_log", "log(%s)");
949 }
950
951 return true;
952 }
953};
954
955class RooExponentialStreamer : public RooFit::JSONIO::Exporter {
956public:
957 std::string const &key() const override;
958 bool exportObject(RooJSONFactoryWSTool *tool, const RooAbsArg *func, JSONNode &elem) const override
959 {
960 auto *pdf = static_cast<const RooExponential *>(func);
961 elem["type"] << key();
962 elem["x"] << pdf->variable().GetName();
963 auto &c = pdf->coefficient();
964 if (pdf->negateCoefficient()) {
965 elem["c"] << c.GetName();
966 } else {
967 elem["c"] << tool->exportTransformed(&c, "_exponential_inverted", "-%s");
968 }
969
970 return true;
971 }
972};
973
974class RooMultiVarGaussianStreamer : public RooFit::JSONIO::Exporter {
975public:
976 std::string const &key() const override;
977 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
978 {
979 auto *pdf = static_cast<const RooMultiVarGaussian *>(func);
980 elem["type"] << key();
981 RooJSONFactoryWSTool::fillSeq(elem["x"], pdf->xVec());
982 RooJSONFactoryWSTool::fillSeq(elem["mean"], pdf->muVec());
983 elem["covariances"].fill_mat(pdf->covarianceMatrix());
984 return true;
985 }
986};
987
988class RooTFnBindingStreamer : public RooFit::JSONIO::Exporter {
989public:
990 std::string const &key() const override;
991 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
992 {
993 auto *pdf = static_cast<const RooTFnBinding *>(func);
994 elem["type"] << key();
995
996 TString formula(pdf->function().GetExpFormula());
997 formula.ReplaceAll("x", pdf->observables()[0].GetName());
998 formula.ReplaceAll("y", pdf->observables()[1].GetName());
999 formula.ReplaceAll("z", pdf->observables()[2].GetName());
1000 for (size_t i = 0; i < pdf->parameters().size(); ++i) {
1001 TString pname(TString::Format("[%d]", (int)i));
1002 formula.ReplaceAll(pname, pdf->parameters()[i].GetName());
1003 }
1004 elem["expression"] << formula.Data();
1005 return true;
1006 }
1007};
1008
1009class RooDerivativeStreamer : public RooFit::JSONIO::Exporter {
1010public:
1011 std::string const &key() const override;
1012 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
1013 {
1014 auto *pdf = static_cast<const RooDerivative *>(func);
1015 elem["type"] << key();
1016 elem["x"] << pdf->getX().GetName();
1017 elem["function"] << pdf->getFunc().GetName();
1018 if (!pdf->getNset().empty()) {
1019 RooJSONFactoryWSTool::fillSeq(elem["normalization"], pdf->getNset());
1020 }
1021 elem["order"] << pdf->order();
1022 elem["eps"] << pdf->eps();
1023 return true;
1024 }
1025};
1026
1027class RooRealIntegralStreamer : public RooFit::JSONIO::Exporter {
1028public:
1029 std::string const &key() const override;
1030 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
1031 {
1032 auto *integral = static_cast<const RooRealIntegral *>(func);
1033 elem["type"] << key();
1034 std::string integrand = integral->integrand().GetName();
1035 // elem["integrand"] << RooJSONFactoryWSTool::sanitizeName(integrand);
1036 elem["integrand"] << integrand;
1037 if (integral->intRange()) {
1038 elem["domain"] << integral->intRange();
1039 }
1040 RooJSONFactoryWSTool::fillSeq(elem["variables"], integral->intVars());
1041 if (RooArgSet const *funcNormSet = integral->funcNormSet()) {
1042 RooJSONFactoryWSTool::fillSeq(elem["normalization"], *funcNormSet);
1043 }
1044 return true;
1045 }
1046};
1047
1048class RooFFTConvPdfStreamer : public RooFit::JSONIO::Exporter {
1049public:
1050 std::string const &key() const override;
1051 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
1052 {
1053 auto *pdf = static_cast<const RooFFTConvPdf *>(func);
1054 elem["type"] << key();
1055 if (auto convFunc = pdf->getPdfConvVar()) {
1056 elem["conv_func"] << convFunc->GetName();
1057 }
1058 elem["conv_var"] << pdf->getConvVar().GetName();
1059 elem["pdf1"] << pdf->getPdf1().GetName();
1060 elem["pdf2"] << pdf->getPdf2().GetName();
1061 elem["ipOrder"] << pdf->getInterpolationOrder();
1062 return true;
1063 }
1064};
1065
1066class RooExtendPdfStreamer : public RooFit::JSONIO::Exporter {
1067public:
1068 std::string const &key() const override;
1069 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
1070 {
1071 auto *pdf = static_cast<const RooExtendPdf *>(func);
1072 elem["type"] << key();
1073 if (auto rangeName = pdf->getRangeName()) {
1074 elem["range"] << rangeName->GetName();
1075 }
1076 elem["pdf"] << pdf->pdf().GetName();
1077 elem["norm"] << pdf->getN().GetName();
1078 return true;
1079 }
1080};
1081
1082class ParamHistFuncStreamer : public RooFit::JSONIO::Exporter {
1083public:
1084 std::string const &key() const override;
1085 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, JSONNode &elem) const override
1086 {
1087 auto *pdf = static_cast<const ParamHistFunc *>(func);
1088 elem["type"] << key();
1089 RooJSONFactoryWSTool::fillSeq(elem["variables"], pdf->dataVars());
1090 RooJSONFactoryWSTool::fillSeq(elem["parameters"], pdf->paramList());
1091 writeBinningInfo(pdf, elem);
1092 return true;
1093 }
1094
1095private:
1096 void writeBinningInfo(const ParamHistFunc *pdf, JSONNode &elem) const
1097 {
1098 auto &observablesNode = elem["axes"].set_seq();
1099 // axes have to be ordered to get consistent bin indices
1100 for (auto *var : static_range_cast<RooRealVar *>(pdf->dataVars())) {
1101 std::string name = var->GetName();
1103 JSONNode &obsNode = observablesNode.append_child().set_map();
1104 obsNode["name"] << name;
1105 if (var->getBinning().isUniform()) {
1106 obsNode["min"] << var->getMin();
1107 obsNode["max"] << var->getMax();
1108 obsNode["nbins"] << var->getBins();
1109 } else {
1110 auto &edges = obsNode["edges"];
1111 edges.set_seq();
1112 double val = var->getBinning().binLow(0);
1113 edges.append_child() << val;
1114 for (int i = 0; i < var->getBinning().numBins(); ++i) {
1115 val = var->getBinning().binHigh(i);
1116 edges.append_child() << val;
1117 }
1118 }
1119 }
1120 }
1121};
1122
1123class RooSplineStreamer : public RooFit::JSONIO::Exporter {
1124public:
1125 std::string const &key() const override;
1126
1127 bool exportObject(RooJSONFactoryWSTool *, const RooAbsArg *func, RooFit::Detail::JSONNode &elem) const override
1128 {
1129 auto const *rs = static_cast<RooSpline const *>(func);
1130
1131 elem["type"] << key();
1132
1133 // Independent variable
1134 elem["x"] << rs->x().GetName();
1135
1136 // Spline configuration
1137 // Canonical algo for RooSpline
1138 elem["interpolation"] << (rs->order() == 5 ? "poly5" : "poly3");
1139 elem["logx"] << rs->logx();
1140 elem["logy"] << rs->logy();
1141
1142 // Serialize knots as primitive arrays
1143 TSpline const &sp = rs->spline();
1144 auto &x0 = elem["x0"].set_seq();
1145 auto &y0 = elem["y0"].set_seq();
1146
1147 const int np = sp.GetNp();
1148 for (int i = 0; i < np; ++i) {
1149 double xk = 0.0, yk = 0.0;
1150 sp.GetKnot(i, xk, yk);
1151 x0.append_child() << xk;
1152 y0.append_child() << yk;
1153 }
1154
1155 return true;
1156 }
1157};
1158
1159#define DEFINE_EXPORTER_KEY(class_name, name) \
1160 std::string const &class_name::key() const \
1161 { \
1162 const static std::string keystring = name; \
1163 return keystring; \
1164 }
1165template <>
1167template <>
1169DEFINE_EXPORTER_KEY(RooBinSamplingPdfStreamer, "binsampling");
1170DEFINE_EXPORTER_KEY(RooBinWidthFunctionStreamer, "binwidth");
1171DEFINE_EXPORTER_KEY(RooLegacyExpPolyStreamer, "legacy_exp_poly_dist");
1172DEFINE_EXPORTER_KEY(RooExponentialStreamer, "exponential_dist");
1173template <>
1175template <>
1177DEFINE_EXPORTER_KEY(RooHistFuncStreamer, "histogram");
1178DEFINE_EXPORTER_KEY(RooHistPdfStreamer, "histogram_dist");
1179DEFINE_EXPORTER_KEY(RooLogNormalStreamer, "lognormal_dist");
1180DEFINE_EXPORTER_KEY(RooMultiVarGaussianStreamer, "multivariate_normal_dist");
1181DEFINE_EXPORTER_KEY(RooPoissonStreamer, "poisson_dist");
1182DEFINE_EXPORTER_KEY(RooDecayStreamer, "decay_dist");
1183DEFINE_EXPORTER_KEY(RooTruthModelStreamer, "truth_model_function");
1184DEFINE_EXPORTER_KEY(RooGaussModelStreamer, "gauss_model_function");
1185template <>
1187template <>
1189DEFINE_EXPORTER_KEY(RooRealSumFuncStreamer, "weighted_sum");
1190DEFINE_EXPORTER_KEY(RooRealSumPdfStreamer, "weighted_sum_dist");
1191DEFINE_EXPORTER_KEY(RooTFnBindingStreamer, "generic_function");
1192DEFINE_EXPORTER_KEY(RooRealIntegralStreamer, "integral");
1193DEFINE_EXPORTER_KEY(RooDerivativeStreamer, "derivative");
1194DEFINE_EXPORTER_KEY(RooFFTConvPdfStreamer, "fft_conv_pdf");
1195DEFINE_EXPORTER_KEY(RooExtendPdfStreamer, "extend_pdf");
1196DEFINE_EXPORTER_KEY(ParamHistFuncStreamer, "step");
1197DEFINE_EXPORTER_KEY(RooSplineStreamer, "spline");
1198
1199///////////////////////////////////////////////////////////////////////////////////////////////////////
1200// instantiate all importers and exporters
1201///////////////////////////////////////////////////////////////////////////////////////////////////////
1202
1203STATIC_EXECUTE([]() {
1204 using namespace RooFit::JSONIO;
1205
1206 registerImporter<RooAddPdfFactory>("mixture_dist", false);
1207 registerImporter<RooAddModelFactory>("mixture_model", false);
1208 registerImporter<RooBinSamplingPdfFactory>("binsampling_dist", false);
1210 registerImporter<RooLegacyExpPolyFactory>("legacy_exp_poly_dist", false);
1211 registerImporter<RooExponentialFactory>("exponential_dist", false);
1212 registerImporter<RooFormulaArgFactory<RooFormulaVar>>("generic_function", false);
1214 registerImporter<RooHistFuncFactory>("histogram", false);
1215 registerImporter<RooHistPdfFactory>("histogram_dist", false);
1216 registerImporter<RooLogNormalFactory>("lognormal_dist", false);
1217 registerImporter<RooMultiVarGaussianFactory>("multivariate_normal_dist", false);
1218 registerImporter<RooPoissonFactory>("poisson_dist", false);
1219 registerImporter<RooDecayFactory>("decay_dist", false);
1220 registerImporter<RooTruthModelFactory>("truth_model_function", false);
1221 registerImporter<RooGaussModelFactory>("gauss_model_function", false);
1224 registerImporter<RooRealSumPdfFactory>("weighted_sum_dist", false);
1225 registerImporter<RooRealSumFuncFactory>("weighted_sum", false);
1227 registerImporter<RooDerivativeFactory>("derivative", false);
1228 registerImporter<RooFFTConvPdfFactory>("fft_conv_pdf", false);
1229 registerImporter<RooExtendPdfFactory>("extend_pdf", false);
1231 registerImporter<RooSplineFactory>("spline", false);
1232
1260});
1261
1262} // namespace
#define DEFINE_EXPORTER_KEY(class_name, name)
bool endsWith(std::string_view str, std::string_view suffix)
std::string removeSuffix(std::string_view str, std::string_view suffix)
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
ROOT::RRangeCast< T, false, Range_t > static_range_cast(Range_t &&coll)
size_t size(const MatrixT &matrix)
retrieve the size of a square matrix
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 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
char name[80]
Definition TGX11.cxx:110
A class which maps the current values of a RooRealVar (or a set of RooRealVars) to one of a number of...
static TClass * Class()
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
bool dependsOn(const RooAbsCollection &serverList, const RooAbsArg *ignoreArg=nullptr, bool valueOnly=false) const
Test whether we depend on (ie, are served by) any object in the specified collection.
virtual bool add(const RooAbsArg &var, bool silent=false)
Add the specified argument to list.
virtual bool addOwned(RooAbsArg &var, bool silent=false)
Add an argument and transfer the ownership to the collection.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
@ CanNotBeExtended
Definition RooAbsPdf.h:208
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
RooAddModel is an efficient implementation of a sum of PDFs of the form.
Definition RooAddModel.h:27
static TClass * Class()
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
static TClass * Class()
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
The RooBinSamplingPdf is supposed to be used as an adapter between a continuous PDF and a binned dist...
static TClass * Class()
double epsilon() const
const RooAbsPdf & pdf() const
const RooAbsReal & observable() const
Returns the bin width (or volume) given a RooHistFunc.
const RooHistFunc & histFunc() const
static TClass * Class()
Implements a RooAbsBinning in terms of an array of boundary values, posing no constraints on the choi...
Definition RooBinning.h:27
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
Single or double sided decay function that can be analytically convolved with any RooResolutionModel ...
Definition RooDecay.h:22
static TClass * Class()
Represents the first, second, or third order derivative of any RooAbsReal as calculated (numerically)...
static TClass * Class()
Exponential PDF.
static TClass * Class()
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
static TClass * Class()
PDF for the numerical (FFT) convolution of two PDFs.
static TClass * Class()
static TClass * Class()
Class RooGaussModel implements a RooResolutionModel that models a Gaussian distribution.
static TClass * Class()
static TClass * Class()
A real-valued function sampled from a multidimensional histogram.
Definition RooHistFunc.h:31
static TClass * Class()
A probability density function sampled from a multidimensional histogram.
Definition RooHistPdf.h:30
static TClass * Class()
When using RooFit, statistical models can be conveniently handled and stored as a RooWorkspace.
static void fillSeq(RooFit::Detail::JSONNode &node, RooAbsCollection const &coll, size_t nMax=-1)
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 bool testValidName(const std::string &str, bool forcError)
static void error(const char *s)
Writes an error message to the RooFit message service and throws a runtime_error.
static std::string name(const RooFit::Detail::JSONNode &n)
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.
RooLegacyExpPoly implements a polynomial PDF of the form.
static TClass * Class()
RooFit Lognormal PDF.
static TClass * Class()
Multivariate Gaussian p.d.f.
static TClass * Class()
Holds the configuration parameters of the various numeric integrators used by RooRealIntegral.
Poisson pdf.
Definition RooPoisson.h:19
static TClass * Class()
static TClass * Class()
static TClass * Class()
Performs hybrid numerical/analytical integrals of RooAbsReal objects.
static TClass * Class()
const RooArgList & coefList() const
const RooArgList & funcList() const
static TClass * Class()
Implements a PDF constructed from a sum of functions:
const RooArgList & funcList() const
static TClass * Class()
ExtendMode extendMode() const override
Returns ability of PDF to provide extended likelihood terms.
const RooArgList & coefList() const
Variable that can be changed from the outside.
Definition RooRealVar.h:37
RooResolutionModel is the base class for PDFs that represent a resolution model that can be convolute...
A RooFit class for creating spline functions.
Definition RooSpline.h:27
static TClass * Class()
Use TF1, TF2, TF3 functions as RooFit objects.
static TClass * Class()
Implements a RooResolution model that corresponds to a delta function.
static TClass * Class()
static Bool_t IsScientificNotation(const TString &formula, int ipos)
Definition TFormula.cxx:345
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
Base class for spline implementation containing the Draw/Paint methods.
Definition TSpline.h:31
Basic string class.
Definition TString.h:138
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
const Double_t sigma
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
#define STATIC_EXECUTE(MY_FUNC)