Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
CodegenImpl.cxx
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 * Authors:
4 * Garima Singh, CERN 2023
5 * Jonas Rembser, CERN 2024
6 *
7 * Copyright (c) 2024, CERN
8 *
9 * Redistribution and use in source and binary forms,
10 * with or without modification, are permitted according to the terms
11 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
12 */
13
14#include <RooFit/CodegenImpl.h>
15
17
18#include <RooAddPdf.h>
19#include <RooAddition.h>
20#include <RooBernstein.h>
21#include <RooBifurGauss.h>
22#include <RooCBShape.h>
23#include <RooCategory.h>
24#include <RooChebychev.h>
25#include <RooConstVar.h>
26#include <RooConstraintSum.h>
27#include <RooEffProd.h>
28#include <RooEfficiency.h>
29#include <RooExponential.h>
30#include <RooExtendPdf.h>
33#include <RooFormulaVar.h>
34#include <RooFunctor1DBinding.h>
35#include <RooFunctorBinding.h>
36#include <RooGamma.h>
37#include <RooGaussian.h>
38#include <RooGenericPdf.h>
39#include <RooHistFunc.h>
40#include <RooHistPdf.h>
41#include <RooLandau.h>
42#include <RooLognormal.h>
43#include <RooMultiPdf.h>
44#include <RooMultiVarGaussian.h>
45#include <RooONNXFunction.h>
46#include <RooParamHistFunc.h>
47#include <RooPoisson.h>
48#include <RooPolyVar.h>
49#include <RooPolynomial.h>
50#include <RooProdPdf.h>
51#include <RooProduct.h>
52#include <RooRatio.h>
53#include <RooRealIntegral.h>
54#include <RooRealSumFunc.h>
55#include <RooRealSumPdf.h>
56#include <RooRealVar.h>
61#include <RooUniform.h>
62#include <RooWrapperPdf.h>
63
64#include "RooFitImplHelpers.h"
65
66#include <TInterpreter.h>
67
68#include <unordered_set>
69
70namespace RooFit::Experimental {
71
72namespace {
73
74// Return a stringy-field version of the value, formatted to maximum precision.
75std::string doubleToString(double val)
76{
77 std::stringstream ss;
78 ss << std::setprecision(std::numeric_limits<double>::max_digits10) << val;
79 return ss.str();
80}
81
82std::string mathFunc(std::string const &name)
83{
84 return "RooFit::Detail::MathFuncs::" + name;
85}
86
87void rooHistTranslateImpl(RooAbsArg const &arg, CodegenContext &ctx, int intOrder, RooDataHist const &dataHist,
88 const RooArgSet &obs, bool correctForBinSize, bool cdfBoundaries)
89{
90 if (intOrder != 0 && !(!cdfBoundaries && !correctForBinSize && intOrder == 1 && obs.size() == 1)) {
91 ooccoutE(&arg, InputArguments) << "RooHistPdf::weight(" << arg.GetName()
92 << ") ERROR: codegen currently only supports non-interpolation cases."
93 << std::endl;
94 return;
95 }
96
97 if (intOrder == 1) {
98 RooAbsBinning const &binning = *dataHist.getBinnings()[0];
99 std::string weightArr = dataHist.declWeightArrayForCodeSquash(ctx, correctForBinSize);
100 ctx.addResult(&arg, ctx.buildCall(mathFunc("interpolate1d"), binning.lowBound(), binning.highBound(), *obs[0],
101 binning.numBins(), weightArr));
102 return;
103 }
104 std::string const &offset = dataHist.calculateTreeIndexForCodeSquash(ctx, obs);
105 std::string weightArr = dataHist.declWeightArrayForCodeSquash(ctx, correctForBinSize);
106 ctx.addResult(&arg, "(" + weightArr + ")[" + offset + "]");
107}
108
109std::string realSumPdfTranslateImpl(CodegenContext &ctx, RooAbsArg const &arg, RooArgList const &funcList,
110 RooArgList const &coefList, bool normalize)
111{
112 bool noLastCoeff = funcList.size() != coefList.size();
113
114 std::string const &funcName = ctx.buildArg(funcList);
115 std::string const &coeffName = ctx.buildArg(coefList);
116 std::string const &coeffSize = std::to_string(coefList.size());
117
118 std::string sum = ctx.getTmpVarName();
119 std::string coeffSum = ctx.getTmpVarName();
120 ctx.addToCodeBody(&arg, "double " + sum + " = 0;\ndouble " + coeffSum + "= 0;\n");
121
122 std::string iterator = "i_" + ctx.getTmpVarName();
123 std::string subscriptExpr = "[" + iterator + "]";
124
125 std::string code = "for(int " + iterator + " = 0; " + iterator + " < " + coeffSize + "; " + iterator + "++) {\n" +
126 sum + " += " + funcName + subscriptExpr + " * " + coeffName + subscriptExpr + ";\n";
127 code += coeffSum + " += " + coeffName + subscriptExpr + ";\n";
128 code += "}\n";
129
130 if (noLastCoeff) {
131 code += sum + " += " + funcName + "[" + coeffSize + "]" + " * (1 - " + coeffSum + ");\n";
132 } else if (normalize) {
133 code += sum + " /= " + coeffSum + ";\n";
134 }
135 ctx.addToCodeBody(&arg, code);
136
137 return sum;
138}
139
140} // namespace
141
143{
144 if (arg.isRearranged()) {
145 ctx.addResult(&arg, ctx.buildCall(mathFunc("ratio"), *arg.rearrangedNum(), *arg.rearrangedDen()));
146 } else {
147 ctx.addResult(&arg, ctx.buildCall(mathFunc("product"), *arg.partList(), arg.partList()->size()));
148 }
149}
150
152{
153 std::string const &idx = arg.dataHist().calculateTreeIndexForCodeSquash(ctx, arg.dataVars(), true);
154 std::string const &paramNames = ctx.buildArg(arg.paramList());
155
156 ctx.addResult(&arg, paramNames + "[" + idx + "]");
157}
158
160{
161 auto const &interpCodes = arg.interpolationCodes();
162
163 std::size_t n = interpCodes.size();
164
165 std::string resName = "total_" + ctx.getTmpVarName();
166 for (std::size_t i = 0; i < n; ++i) {
167 if (interpCodes[i] != interpCodes[0]) {
169 << "FlexibleInterpVar::evaluate ERROR: Code Squashing AD does not yet support having "
170 "different interpolation codes for the same class object "
171 << std::endl;
172 }
173 }
174
175 // The PiecewiseInterpolation class is used in the context of HistFactory
176 // models, where is is always used the same way: all RooAbsReals in _lowSet,
177 // _histSet, and also nominal are 1D RooHistFuncs with with same structure.
178 //
179 // Therefore, we can make a big optimization: we get the bin index only once
180 // here in the generated code for PiecewiseInterpolation. Then, we also
181 // rearrange the histogram data in such a way that we can always pass the
182 // same arrays to the free function that implements the interpolation, just
183 // with a dynamic offset calculated from the bin index.
184 RooDataHist const &nomHist = dynamic_cast<RooHistFunc const &>(*arg.nominalHist()).dataHist();
185 int nBins = nomHist.numEntries();
186 std::vector<double> valsNominal;
187 std::vector<double> valsLow;
188 std::vector<double> valsHigh;
189 for (int i = 0; i < nBins; ++i) {
190 valsNominal.push_back(nomHist.weight(i));
191 }
192 for (int i = 0; i < nBins; ++i) {
193 for (std::size_t iParam = 0; iParam < n; ++iParam) {
194 valsLow.push_back(dynamic_cast<RooHistFunc const &>(arg.lowList()[iParam]).dataHist().weight(i));
195 valsHigh.push_back(dynamic_cast<RooHistFunc const &>(arg.highList()[iParam]).dataHist().weight(i));
196 }
197 }
198 std::string idxName = ctx.getTmpVarName();
199 std::string valsNominalStr = ctx.buildArg(valsNominal);
200 std::string valsLowStr = ctx.buildArg(valsLow);
201 std::string valsHighStr = ctx.buildArg(valsHigh);
202 std::string nStr = std::to_string(n);
203 std::string code;
204
205 std::string lowName = ctx.getTmpVarName();
206 std::string highName = ctx.getTmpVarName();
207 std::string nominalName = ctx.getTmpVarName();
208 code +=
209 "unsigned int " + idxName + " = " +
210 nomHist.calculateTreeIndexForCodeSquash(ctx, dynamic_cast<RooHistFunc const &>(*arg.nominalHist()).variables()) +
211 ";\n";
212 code += "double const* " + lowName + " = " + valsLowStr + " + " + nStr + " * " + idxName + ";\n";
213 code += "double const* " + highName + " = " + valsHighStr + " + " + nStr + " * " + idxName + ";\n";
214 code += "double " + nominalName + " = *(" + valsNominalStr + " + " + idxName + ");\n";
215
216 std::string funcCall = ctx.buildCall(mathFunc("flexibleInterp"), interpCodes[0], arg.paramList(), n, lowName,
217 highName, 1.0, nominalName, 0.0);
218 code += "double " + resName + " = " + funcCall + ";\n";
219
220 if (arg.positiveDefinite()) {
221 code += resName + " = " + resName + " < 0 ? 0 : " + resName + ";\n";
222 }
223
224 ctx.addToCodeBody(&arg, code);
225 ctx.addResult(&arg, resName);
226}
227
228////////////////////////////////////////////////////////////////////////////////
229/// This function defines a translation for each RooAbsReal based object that can be used
230/// to express the class as simple C++ code. The function adds the code represented by
231/// each class as an std::string (that is later concatenated with code strings from translate calls)
232/// to form the C++ code that AD tools can understand. Any class that wants to support AD, has to
233/// implement this function.
234///
235/// \param[in] ctx An object to manage auxiliary information for code-squashing. Also takes the
236/// code string that this class outputs into the squashed code through the 'addToCodeBody' function.
238{
239 std::stringstream errorMsg;
240 errorMsg << "Translate function for class \"" << arg.ClassName() << "\" has not yet been implemented.";
241 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
242 return ctx.addResult(&arg, "1.0");
243}
244
246{
247 ctx.addResult(&arg, realSumPdfTranslateImpl(ctx, arg, arg.pdfList(), arg.coefList(), true));
248}
249
251{
252 auto const &covI = arg.covarianceMatrixInverse();
253 std::span<const double> covISpan{covI.GetMatrixArray(), static_cast<size_t>(covI.GetNoElements())};
254 ctx.addResult(&arg,
255 ctx.buildCall(mathFunc("multiVarGaussian"), arg.xVec().size(), arg.xVec(), arg.muVec(), covISpan));
256}
257
259{
260 int numPdfs = arg.getNumPdfs();
261
262 // MathFunc call
263
264 // The value of this number should be discussed. Beyound a certain number of
265 // indices MathFunc call becomes more efficient.
266 if (numPdfs > 2) {
267 ctx.addResult(&arg, ctx.buildCall(mathFunc("multipdf"), arg.indexCategory(), arg.getPdfList()));
268
269 std::cout << "MathFunc call used\n";
270
271 } else {
272
273 // Ternary nested expression
274 std::string indexExpr = ctx.getResult(arg.indexCategory());
275
276 // int numPdfs = arg.getNumPdfs();
277 std::string expr;
278
279 for (int i = 0; i < numPdfs; ++i) {
280 RooAbsPdf *pdf = arg.getPdf(i);
281 std::string pdfExpr = ctx.getResult(*pdf);
282
283 expr += "(" + indexExpr + " == " + std::to_string(i) + " ? (" + pdfExpr + ") : ";
284 }
285
286 expr += "0.0";
287 expr += std::string(numPdfs, ')'); // Close all ternary operators
288
289 ctx.addResult(&arg, expr);
290 std::cout << "Ternary expression call used \n";
291 }
292}
293
294// RooCategory index added.
296{
297 int idx = ctx.observableIndexOf(arg);
298 if (idx < 0) {
299
300 idx = 1;
301 ctx.addVecObs(arg.GetName(), idx);
302 }
303
304 std::string result = std::to_string(arg.getCurrentIndex());
305 ctx.addResult(&arg, result);
306}
307
309{
310 if (arg.list().empty()) {
311 ctx.addResult(&arg, "0.0");
312 }
313 std::string result;
314 if (arg.list().size() > 1)
315 result += "(";
316
317 std::size_t i = 0;
318 for (auto *component : static_range_cast<RooAbsReal *>(arg.list())) {
319
320 if (!dynamic_cast<RooFit::Detail::RooNLLVarNew *>(component) || arg.list().size() == 1) {
321 result += ctx.getResult(*component);
322 ++i;
323 if (i < arg.list().size())
324 result += '+';
325 continue;
326 }
327 result += ctx.buildFunction(*component, ctx.dependsOnData()) + "(params, obs, xlArr)";
328 ++i;
329 if (i < arg.list().size())
330 result += '+';
331 }
332 if (arg.list().size() > 1)
333 result += ')';
334 ctx.addResult(&arg, result);
335}
336
338{
339 arg.fillBuffer();
340 ctx.addResult(&arg, ctx.buildCall(mathFunc("bernstein"), arg.x(), arg.xmin(), arg.xmax(), arg.coefList(),
341 arg.coefList().size()));
342}
343
345{
346 ctx.addResult(&arg,
347 ctx.buildCall(mathFunc("bifurGauss"), arg.getX(), arg.getMean(), arg.getSigmaL(), arg.getSigmaR()));
348}
349
351{
352 ctx.addResult(
353 &arg, ctx.buildCall(mathFunc("cbShape"), arg.getM(), arg.getM0(), arg.getSigma(), arg.getAlpha(), arg.getN()));
354}
355
357{
358 // first bring the range of the variable _x to the normalised range [-1, 1]
359 // calculate sum_k c_k T_k(x) where x is given in the normalised range,
360 // c_0 = 1, and the higher coefficients are given in _coefList
361 double xmax = static_cast<RooAbsRealLValue const &>(arg.x()).getMax(arg.refRangeName());
362 double xmin = static_cast<RooAbsRealLValue const &>(arg.x()).getMin(arg.refRangeName());
363
364 ctx.addResult(&arg,
365 ctx.buildCall(mathFunc("chebychev"), arg.coefList(), arg.coefList().size(), arg.x(), xmin, xmax));
366}
367
369{
370 ctx.addResult(&arg, doubleToString(arg.getVal()));
371}
372
374{
375 ctx.addResult(&arg, ctx.buildCall(mathFunc("constraintSum"), arg.list(), arg.list().size()));
376}
377
378// Generate RooFit codegen wrappers for RooFunctorBinding and similar objects,
379// emitting both the primal function call and its gradient pullback for
380// Clad-based AD.
381template <class RooArg_t>
382void functorCodegenImpl(RooArg_t &arg, RooArgList const &variables, CodegenContext &ctx)
383{
384 if (!arg.function()->HasGradient()) {
385 std::stringstream errorMsg;
386 errorMsg << "Functor wrapped by \"" << arg.GetName() << "\" doesn't provide a gradient function."
387 << " RooFit codegen is therefore not supported.";
388 oocoutE(&arg, InputArguments) << errorMsg.str() << std::endl;
389 throw std::runtime_error(errorMsg.str());
390 }
391
392 std::string funcAddrStr = TString::Format("0x%zx", reinterpret_cast<std::size_t>(arg.function())).Data();
393 std::string wrapperName = "roo_functor_" + funcAddrStr;
394
395 static std::unordered_set<std::string> wrapperNames;
396
397 if (wrapperNames.find(wrapperName) == wrapperNames.end()) {
398
400
401 std::string pullbackName = wrapperName + "_pullback";
402 std::string nStr = std::to_string(std::size(variables));
403
404 std::string type;
405 if constexpr (std::is_same_v<RooArg_t, RooFunctor1DBinding> || std::is_same_v<RooArg_t, RooFunctor1DPdfBinding>)
406 type = "::ROOT::Math::IGradientFunctionOneDim";
407 else
408 type = "::ROOT::Math::IGradientFunctionMultiDim";
409
410 std::string funcAddrCasted = "reinterpret_cast<" + type + " const *>(" + funcAddrStr + ")";
411
412 std::string code;
413
414 code += "double " + wrapperName +
415 "(double const *x) {\n"
416 " return " +
418 "->operator()(x);\n"
419 "}\n\n"
420 "namespace clad::custom_derivatives {\n\n"
421 "void " +
423 "(double const* x, double d_y, double *d_x) {\n"
424 " double output[" +
425 nStr +
426 "]{};\n"
427 " " +
429 "->Gradient(x, output);\n"
430 " for (int i = 0; i < " +
431 nStr +
432 "; ++i) {\n"
433 " d_x[i] += output[i] * d_y;\n"
434 " }\n"
435 "}\n"
436 "} // namespace clad::custom_derivatives\n";
437
438 gInterpreter->Declare(code.c_str());
439 }
440
441 ctx.addResult(&arg, ctx.buildCall(wrapperName, variables));
442}
443
445{
446 functorCodegenImpl(arg, arg.variable(), ctx);
447}
448
450{
451 functorCodegenImpl(arg, arg.variable(), ctx);
452}
453
455{
456 functorCodegenImpl(arg, arg.variables(), ctx);
457}
458
460{
461 functorCodegenImpl(arg, arg.variables(), ctx);
462}
463
465{
466 ctx.addResult(&arg, ctx.buildCall("TMath::GammaDist", arg.getX(), arg.getGamma(), arg.getMu(), arg.getBeta()));
467}
468
470{
471 arg.getVal(); // to trigger the creation of the TFormula
472 std::string funcName = arg.getUniqueFuncName();
473 ctx.collectFunction(funcName);
474 // We have to force the array type to be "double" because that's what the
475 // declared function wrapped by the TFormula expects.
476 auto inputVar = ctx.buildArg(arg.dependents(), /*arrayType=*/"double");
477 ctx.addResult(&arg, funcName + "(" + inputVar + ")");
478}
479
481{
482 ctx.addResult(&arg, ctx.buildCall(mathFunc("effProd"), arg.eff(), arg.pdf()));
483}
484
486{
487 RooAbsCategory const &cat = arg.cat();
488 int sigCatIndex = cat.lookupIndex(arg.sigCatName());
489 ctx.addResult(&arg, ctx.buildCall(mathFunc("efficiency"), arg.effFunc(), cat, sigCatIndex));
490}
491
493{
494 // Build a call to the stateless exponential defined later.
495 std::string coef;
496 if (arg.negateCoefficient()) {
497 coef += "-";
498 }
499 coef += ctx.getResult(arg.coefficient());
500 ctx.addResult(&arg, "std::exp(" + coef + " * " + ctx.getResult(arg.variable()) + ")");
501}
502
504{
505 // Use the result of the underlying pdf.
506 ctx.addResult(&arg, ctx.getResult(arg.pdf()));
507}
508
510{
511 // Build a call to the stateless gaussian defined later.
512 ctx.addResult(&arg, ctx.buildCall(mathFunc("gaussian"), arg.getX(), arg.getMean(), arg.getSigma()));
513}
514
516{
517 arg.getVal(); // to trigger the creation of the TFormula
518 std::string funcName = arg.getUniqueFuncName();
519 ctx.collectFunction(funcName);
520 // We have to force the array type to be "double" because that's what the
521 // declared function wrapped by the TFormula expects.
522 auto inputVar = ctx.buildArg(arg.dependents(), /*arrayType=*/"double");
523 ctx.addResult(&arg, funcName + "(" + inputVar + ")");
524}
525
527{
528 rooHistTranslateImpl(arg, ctx, arg.getInterpolationOrder(), arg.dataHist(), arg.variables(), false,
529 arg.getCdfBoundaries());
530}
531
533{
534 rooHistTranslateImpl(arg, ctx, arg.getInterpolationOrder(), arg.dataHist(), arg.variables(), !arg.haveUnitNorm(),
535 arg.getCdfBoundaries());
536}
537
539{
540 ctx.addResult(&arg, ctx.buildCall(mathFunc("landau"), arg.getX(), arg.getMean(), arg.getSigma()));
541}
542
544{
545 std::string funcName = arg.useStandardParametrization() ? "logNormalEvaluateStandard" : "logNormal";
546 ctx.addResult(&arg, ctx.buildCall(mathFunc(funcName), arg.getX(), arg.getShapeK(), arg.getMedian()));
547}
548
549void codegenImpl(RooFit::Detail::RooNLLVarNew &arg, CodegenContext &ctx)
550{
551 if (arg.binnedL() && !arg.pdf().getAttribute("BinnedLikelihoodActiveYields")) {
552 std::stringstream errorMsg;
553 errorMsg << "codegen: binned likelihood optimization is only supported when raw pdf "
554 "values can be interpreted as yields."
555 << " This is not the case for HistFactory models written with ROOT versions before 6.26.00";
556 oocoutE(&arg, InputArguments) << errorMsg.str() << std::endl;
557 throw std::runtime_error(errorMsg.str());
558 }
559
560 std::string weightSumName = RooFit::Detail::makeValidVarName(arg.GetName()) + "WeightSum";
561 std::string resName = RooFit::Detail::makeValidVarName(arg.GetName()) + "Result";
562 ctx.addResult(&arg, resName);
563 ctx.addToGlobalScope("double " + weightSumName + " = 0.0;\n");
564 ctx.addToGlobalScope("double " + resName + " = 0.0;\n");
565
566 const bool needWeightSum = arg.expectedEvents() || arg.simCount() > 1;
567
568 if (needWeightSum) {
569 auto scope = ctx.beginLoop(&arg);
570 ctx.addToCodeBody(weightSumName + " += " + ctx.getResult(arg.weightVar()) + ";\n");
571 }
572 if (arg.simCount() > 1) {
573 std::string simCountStr = std::to_string(static_cast<double>(arg.simCount()));
574 ctx.addToCodeBody(resName + " += " + weightSumName + " * std::log(" + simCountStr + ");\n");
575 }
576
577 // Begin loop scope for the observables and weight variable. If the weight
578 // is a scalar, the context will ignore it for the loop scope. The closing
579 // brackets of the loop is written at the end of the scopes lifetime.
580 {
581 auto scope = ctx.beginLoop(&arg);
582 std::string term = ctx.buildCall(mathFunc("nll"), arg.pdf(), arg.weightVar(), arg.binnedL(), 0);
583 ctx.addToCodeBody(&arg, resName + " += " + term + ";");
584 }
585 if (arg.expectedEvents()) {
586 std::string expected = ctx.getResult(*arg.expectedEvents());
587 ctx.addToCodeBody(resName + " += " + expected + " - " + weightSumName + " * std::log(" + expected + ");\n");
588 }
589}
590
592{
593 // For now just return function/normalization integral.
594 ctx.addResult(&arg, ctx.getResult(arg.pdf()) + "/" + ctx.getResult(arg.normIntegral()));
595}
596
598{
599 std::string const &idx = arg.dataHist().calculateTreeIndexForCodeSquash(ctx, arg.xList());
600 std::string arrName = ctx.buildArg(arg.paramList());
601 std::stringstream result;
602 result << arrName << "[" << idx << "]";
603 if (arg.relParam()) {
604 // get weight[idx] * binv[idx]. Here we get the bin volume for the first element as we assume the distribution to
605 // be binned uniformly.
606 double binV = arg.dataHist().binVolume(0);
607 std::string weightArr = arg.dataHist().declWeightArrayForCodeSquash(ctx, false);
608 result << " * *(" << weightArr << " + " << idx + ") * " << doubleToString(binV);
609 }
610 ctx.addResult(&arg, result.str());
611}
612
614{
615 std::string xName = ctx.getResult(arg.getX());
616 if (!arg.getNoRounding())
617 xName = "std::floor(" + xName + ")";
618
619 ctx.addResult(&arg, ctx.buildCall(mathFunc("poisson"), xName, arg.getMean()));
620}
621
623{
624 const unsigned sz = arg.coefList().size();
625 if (!sz) {
626 ctx.addResult(&arg, std::to_string(arg.lowestOrder() ? 1. : 0.));
627 return;
628 }
629
630 ctx.addResult(&arg, ctx.buildCall(mathFunc("polynomial"), arg.coefList(), sz, arg.lowestOrder(), arg.x()));
631}
632
634{
635 const unsigned sz = arg.coefList().size();
636 if (!sz) {
637 ctx.addResult(&arg, std::to_string(arg.lowestOrder() ? 1. : 0.));
638 return;
639 }
640
641 ctx.addResult(&arg, ctx.buildCall(mathFunc("polynomial<true>"), arg.coefList(), sz, arg.lowestOrder(), arg.x()));
642}
643
645{
646 ctx.addResult(&arg, ctx.buildCall(mathFunc("product"), arg.realComponents(), arg.realComponents().size()));
647}
648
650{
651 ctx.addResult(&arg, ctx.buildCall(mathFunc("ratio"), arg.numerator(), arg.denominator()));
652}
653
654namespace {
655
656std::string codegenIntegral(RooAbsReal &arg, int code, const char *rangeName, CodegenContext &ctx)
657{
658 using Func = std::string (*)(RooAbsReal &, int, const char *, CodegenContext &);
659
660 Func func;
661
662 TClass *tclass = arg.IsA();
663
664 // Cache the overload resolutions
665 static std::unordered_map<TClass *, Func> dispatchMap;
666
667 auto found = dispatchMap.find(tclass);
668
669 if (found != dispatchMap.end()) {
670 func = found->second;
671 } else {
672 // Can probably done with CppInterop in the future to avoid string manipulation.
673 std::stringstream cmd;
674 cmd << "&RooFit::Experimental::CodegenIntegralImplCaller<" << tclass->GetName() << ">::call;";
675 func = reinterpret_cast<Func>(gInterpreter->ProcessLine(cmd.str().c_str()));
676 dispatchMap[tclass] = func;
677 }
678
679 return func(arg, code, rangeName, ctx);
680}
681
682} // namespace
683
685{
686 if (arg.numIntCatVars().empty() && arg.numIntRealVars().empty()) {
687 ctx.addResult(&arg, codegenIntegral(const_cast<RooAbsReal &>(arg.integrand()), arg.mode(), arg.intRange(), ctx));
688 return;
689 }
690
691 if (arg.intVars().size() != 1 || arg.numIntRealVars().size() != 1) {
692 std::stringstream errorMsg;
693 errorMsg << "Only analytical integrals and 1D numeric integrals are supported for AD for class"
694 << arg.integrand().GetName();
695 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
696 throw std::runtime_error(errorMsg.str().c_str());
697 }
698
699 auto &intVar = static_cast<RooAbsRealLValue &>(*arg.numIntRealVars()[0]);
700
701 std::string obsName = ctx.getTmpVarName();
702 std::string oldIntVarResult = ctx.getResult(intVar);
703 ctx.addResult(&intVar, "obs[0]");
704
705 std::string funcName = ctx.buildFunction(arg.integrand(), {});
706
707 std::stringstream ss;
708
709 ss << "double " << obsName << "[1];\n";
710
711 std::string resName = RooFit::Detail::makeValidVarName(arg.GetName()) + "Result";
712 ctx.addResult(&arg, resName);
713 ctx.addToGlobalScope("double " + resName + " = 0.0;\n");
714
715 // TODO: once Clad has support for higher-order functions (follow also the
716 // Clad issue #637), we could refactor this code into an actual function
717 // instead of hardcoding it here as a string.
718 ss << "{\n"
719 << " const int n = 1000; // number of sampling points\n"
720 << " double d = " << intVar.getMax(arg.intRange()) << " - " << intVar.getMin(arg.intRange()) << ";\n"
721 << " double eps = d / n;\n"
722 << " #pragma clad checkpoint loop\n"
723 << " for (int i = 0; i < n; ++i) {\n"
724 << " " << obsName << "[0] = " << intVar.getMin(arg.intRange()) << " + eps * i;\n"
725 << " double tmpA = " << funcName << "(params, " << obsName << ", xlArr);\n"
726 << " " << obsName << "[0] = " << intVar.getMin(arg.intRange()) << " + eps * (i + 1);\n"
727 << " double tmpB = " << funcName << "(params, " << obsName << ", xlArr);\n"
728 << " " << resName << " += (tmpA + tmpB) * 0.5 * eps;\n"
729 << " }\n"
730 << "}\n";
731
732 ctx.addToGlobalScope(ss.str());
733
735}
736
738{
739 ctx.addResult(&arg, realSumPdfTranslateImpl(ctx, arg, arg.funcList(), arg.coefList(), false));
740}
741
743{
744 ctx.addResult(&arg, realSumPdfTranslateImpl(ctx, arg, arg.funcList(), arg.coefList(), false));
745}
746
748{
749 if (!arg.isConstant()) {
750 ctx.addResult(&arg, arg.GetName());
751 }
752 ctx.addResult(&arg, doubleToString(arg.getVal()));
753}
754
756{
757 ctx.addResult(&arg, ctx.buildCall(mathFunc("recursiveFraction"), arg.variables(), arg.variables().size()));
758}
759
761{
762 auto const &interpCodes = arg.interpolationCodes();
763
764 unsigned int n = interpCodes.size();
765
766 int interpCode = interpCodes[0];
767 // To get consistent codes with the PiecewiseInterpolation
768 if (interpCode == 4) {
769 interpCode = 5;
770 }
771
772 for (unsigned int i = 1; i < n; i++) {
773 if (interpCodes[i] != interpCodes[0]) {
775 << "FlexibleInterpVar::evaluate ERROR: Code Squashing AD does not yet support having "
776 "different interpolation codes for the same class object "
777 << std::endl;
778 }
779 }
780
781 std::string const &resName = ctx.buildCall(mathFunc("flexibleInterp"), interpCode, arg.variables(), n, arg.low(),
782 arg.high(), arg.globalBoundary(), arg.nominal(), 1.0);
783 ctx.addResult(&arg, resName);
784}
785
787{
788 ctx.addResult(&arg, "1.0");
789}
790
792{
793 ctx.addResult(&arg, ctx.getResult(arg.function()));
794}
795
796////////////////////////////////////////////////////////////////////////////////
797/// This function defines the analytical integral translation for the class.
798///
799/// \param[in] code The code that decides the integrands.
800/// \param[in] rangeName Name of the normalization range.
801/// \param[in] ctx An object to manage auxiliary information for code-squashing.
802///
803/// \returns The representative code string of the integral for the given object.
804std::string codegenIntegralImpl(RooAbsReal &arg, int, const char *, CodegenContext &)
805{
806 std::stringstream errorMsg;
807 errorMsg << "An analytical integral function for class \"" << arg.ClassName() << "\" has not yet been implemented.";
808 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
809 throw std::runtime_error(errorMsg.str().c_str());
810}
811
812std::string codegenIntegralImpl(RooBernstein &arg, int, const char *rangeName, CodegenContext &ctx)
813{
814 arg.fillBuffer(); // to get the right xmin() and xmax()
815 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
816 return ctx.buildCall(mathFunc("bernsteinIntegral"), x.getMin(rangeName), x.getMax(rangeName), arg.xmin(), arg.xmax(),
817 arg.coefList(), arg.coefList().size());
818}
819
820std::string codegenIntegralImpl(RooBifurGauss &arg, int code, const char *rangeName, CodegenContext &ctx)
821{
822 auto &constant = code == 1 ? arg.getMean() : arg.getX();
823 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(code == 1 ? arg.getX() : arg.getMean());
824
825 return ctx.buildCall(mathFunc("bifurGaussIntegral"), integrand.getMin(rangeName), integrand.getMax(rangeName),
826 constant, arg.getSigmaL(), arg.getSigmaR());
827}
828
829std::string codegenIntegralImpl(RooCBShape &arg, int /*code*/, const char *rangeName, CodegenContext &ctx)
830{
831 auto &m = dynamic_cast<RooAbsRealLValue const &>(arg.getM());
832 return ctx.buildCall(mathFunc("cbShapeIntegral"), m.getMin(rangeName), m.getMax(rangeName), arg.getM0(),
833 arg.getSigma(), arg.getAlpha(), arg.getN());
834}
835
836std::string codegenIntegralImpl(RooChebychev &arg, int, const char *rangeName, CodegenContext &ctx)
837{
838 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
839 double xmax = x.getMax(arg.refRangeName());
840 double xmin = x.getMin(arg.refRangeName());
841 unsigned int sz = arg.coefList().size();
842
843 return ctx.buildCall(mathFunc("chebychevIntegral"), arg.coefList(), sz, xmin, xmax, x.getMin(rangeName),
844 x.getMax(rangeName));
845}
846
847std::string codegenIntegralImpl(RooEfficiency &, int, const char *, CodegenContext &)
848{
849 return "1.0";
850}
851
852std::string codegenIntegralImpl(RooExponential &arg, int code, const char *rangeName, CodegenContext &ctx)
853{
854 bool isOverX = code == 1;
855
856 std::string constant;
857 if (arg.negateCoefficient() && isOverX) {
858 constant += "-";
859 }
860 constant += ctx.getResult(isOverX ? arg.coefficient() : arg.variable());
861
862 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(isOverX ? arg.variable() : arg.coefficient());
863
864 double min = integrand.getMin(rangeName);
865 double max = integrand.getMax(rangeName);
866
867 if (!isOverX && arg.negateCoefficient()) {
868 std::swap(min, max);
869 min = -min;
870 max = -max;
871 }
872
873 return ctx.buildCall(mathFunc("exponentialIntegral"), min, max, constant);
874}
875
876std::string codegenIntegralImpl(RooGamma &arg, int, const char *rangeName, CodegenContext &ctx)
877{
878 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.getX());
879 const std::string a =
880 ctx.buildCall("ROOT::Math::gamma_cdf", x.getMax(rangeName), arg.getGamma(), arg.getBeta(), arg.getMu());
881 const std::string b =
882 ctx.buildCall("ROOT::Math::gamma_cdf", x.getMin(rangeName), arg.getGamma(), arg.getBeta(), arg.getMu());
883 return a + " - " + b;
884}
885
886std::string codegenIntegralImpl(RooGaussian &arg, int code, const char *rangeName, CodegenContext &ctx)
887{
888 auto &constant = code == 1 ? arg.getMean() : arg.getX();
889 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(code == 1 ? arg.getX() : arg.getMean());
890
891 return ctx.buildCall(mathFunc("gaussianIntegral"), integrand.getMin(rangeName), integrand.getMax(rangeName),
892 constant, arg.getSigma());
893}
894
895namespace {
896
897std::string rooHistIntegralTranslateImpl(int code, RooAbsArg const &arg, RooDataHist const &dataHist,
898 const RooArgSet &obs, bool histFuncMode)
899{
900 if (((2 << obs.size()) - 1) != code) {
901 oocoutE(&arg, InputArguments) << "RooHistPdf::integral(" << arg.GetName()
902 << ") ERROR: AD currently only supports integrating over all histogram observables."
903 << std::endl;
904 return "";
905 }
906 return doubleToString(dataHist.sum(histFuncMode));
907}
908
909} // namespace
910
911std::string codegenIntegralImpl(RooHistFunc &arg, int code, const char *, CodegenContext &)
912{
913 return rooHistIntegralTranslateImpl(code, arg, arg.dataHist(), arg.variables(), true);
914}
915
916std::string codegenIntegralImpl(RooHistPdf &arg, int code, const char *, CodegenContext &)
917{
918 return rooHistIntegralTranslateImpl(code, arg, arg.dataHist(), arg.variables(), false);
919}
920
921std::string codegenIntegralImpl(RooLandau &arg, int, const char *rangeName, CodegenContext &ctx)
922{
923 // Don't do anything with "code". It can only be "1" anyway (see
924 // implementation of getAnalyticalIntegral).
925 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.getX());
926 const std::string a = ctx.buildCall("ROOT::Math::landau_cdf", x.getMax(rangeName), arg.getSigma(), arg.getMean());
927 const std::string b = ctx.buildCall("ROOT::Math::landau_cdf", x.getMin(rangeName), arg.getSigma(), arg.getMean());
928 return ctx.getResult(arg.getSigma()) + " * " + "(" + a + " - " + b + ")";
929}
930
931std::string codegenIntegralImpl(RooLognormal &arg, int, const char *rangeName, CodegenContext &ctx)
932{
933 std::string funcName = arg.useStandardParametrization() ? "logNormalIntegralStandard" : "logNormalIntegral";
934 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.getX());
935 return ctx.buildCall(mathFunc(funcName), x.getMin(rangeName), x.getMax(rangeName), arg.getMedian(), arg.getShapeK());
936}
937
938std::string codegenIntegralImpl(RooMultiVarGaussian &arg, int code, const char *rangeName, CodegenContext &)
939{
940 if (code != -1) {
941 std::stringstream errorMsg;
942 errorMsg << "Partial integrals over RooMultiVarGaussian are not supported.";
943 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
944 throw std::runtime_error(errorMsg.str().c_str());
945 }
946
948}
949
951{
952 std::stringstream ss;
953 ss << arg.outerWrapperName() << "(";
954 for (std::size_t i = 0; i < arg.nInputTensors(); ++i) {
955 ss << ctx.buildArg(arg.inputTensorList(i)) << std::endl;
956 if (i != arg.nInputTensors() - 1) {
957 ss << ", ";
958 }
959 }
960 ss << ")";
961
962 ctx.addResult(&arg, ss.str());
963}
964
965std::string codegenIntegralImpl(RooPoisson &arg, int code, const char *rangeName, CodegenContext &ctx)
966{
967 assert(code == 1 || code == 2);
968 std::string xName = ctx.getResult(arg.getX());
969 if (!arg.getNoRounding())
970 xName = "std::floor(" + xName + ")";
971
972 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(code == 1 ? arg.getX() : arg.getMean());
973 // Since the integral function is the same for both codes, we need to make sure the indexed observables do not appear
974 // in the function if they are not required.
975 xName = code == 1 ? "0" : xName;
976 return ctx.buildCall(mathFunc("poissonIntegral"), code, arg.getMean(), xName, integrand.getMin(rangeName),
977 integrand.getMax(rangeName), arg.getProtectNegativeMean());
978}
979
980std::string codegenIntegralImpl(RooPolyVar &arg, int, const char *rangeName, CodegenContext &ctx)
981{
982 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
983 const double xmin = x.getMin(rangeName);
984 const double xmax = x.getMax(rangeName);
985 const unsigned sz = arg.coefList().size();
986 if (!sz)
987 return std::to_string(arg.lowestOrder() ? xmax - xmin : 0.0);
988
989 return ctx.buildCall(mathFunc("polynomialIntegral"), arg.coefList(), sz, arg.lowestOrder(), xmin, xmax);
990}
991
992std::string codegenIntegralImpl(RooPolynomial &arg, int, const char *rangeName, CodegenContext &ctx)
993{
994 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
995 const double xmin = x.getMin(rangeName);
996 const double xmax = x.getMax(rangeName);
997 const unsigned sz = arg.coefList().size();
998 if (!sz)
999 return std::to_string(arg.lowestOrder() ? xmax - xmin : 0.0);
1000
1001 return ctx.buildCall(mathFunc("polynomialIntegral<true>"), arg.coefList(), sz, arg.lowestOrder(), xmin, xmax);
1002}
1003
1004std::string codegenIntegralImpl(RooRealSumPdf &arg, int code, const char *rangeName, CodegenContext &ctx)
1005{
1006 // Re-use translate, since integration is linear.
1007 return realSumPdfTranslateImpl(ctx, arg, arg.funcIntListFromCache(code, rangeName), arg.coefList(), false);
1008}
1009
1010std::string codegenIntegralImpl(RooUniform &arg, int code, const char *rangeName, CodegenContext &)
1011{
1012 // The integral of a uniform distribution is static, so we can just hardcode
1013 // the result in a string.
1014 return doubleToString(arg.analyticalIntegral(code, rangeName));
1015}
1016
1017} // namespace RooFit::Experimental
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#define oocoutE(o, a)
#define ooccoutE(o, a)
ROOT::Detail::TRangeCast< T, true > TRangeDynCast
TRangeDynCast is an adapter class that allows the typed iteration through a TCollection.
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 offset
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 TString Int_t GCValues_t GetPrimarySelectionOwner GetDisplay GetScreen GetColormap GetNativeEvent const char const char dpyName wid window const char font_name cursor keysym reg const char only_if_exist regb h Point_t winding char text const char depth char const char Int_t count const char ColorStruct_t color const char Pixmap_t Pixmap_t PictureAttributes_t attr const char char ret_data h unsigned char height h Atom_t Int_t ULong_t ULong_t unsigned char prop_list Atom_t Atom_t Atom_t Time_t type
char name[80]
Definition TGX11.cxx:148
float xmin
float xmax
#define gInterpreter
A class which maps the current values of a RooRealVar (or a set of RooRealVars) to one of a number of...
const RooArgList & paramList() const
const RooArgList & dataVars() const
RooDataHist const & dataHist() const
The PiecewiseInterpolation is a class that can morph distributions into each other,...
const RooArgList & highList() const
const RooAbsReal * nominalHist() const
Return pointer to the nominal hist function.
const RooArgList & lowList() const
const RooArgList & paramList() const
const std::vector< int > & interpolationCodes() 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 isConstant() const
Check if the "Constant" attribute is set.
Definition RooAbsArg.h:283
Abstract base class for RooRealVar binning definitions.
Int_t numBins() const
Return number of bins.
virtual double highBound() const =0
virtual double lowBound() const =0
A space to attach TBranches.
value_type lookupIndex(const std::string &stateName) const
Find the index number corresponding to the state name.
Storage_t::size_type size() const
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:32
Abstract base class for objects that represent a real value that may appear on the left hand side of ...
virtual double getMin(const char *name=nullptr) const
Get minimum of currently defined range.
Abstract base class for objects that represent a real value and implements functionality common to al...
Definition RooAbsReal.h:63
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:107
TClass * IsA() const override
Definition RooAbsReal.h:551
Efficient implementation of a sum of PDFs of the form.
Definition RooAddPdf.h:33
const RooArgList & coefList() const
Definition RooAddPdf.h:74
const RooArgList & pdfList() const
Definition RooAddPdf.h:70
Calculates the sum of a set of RooAbsReal terms, or when constructed with two sets,...
Definition RooAddition.h:27
const RooArgList & list() const
Definition RooAddition.h:42
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
Bernstein basis polynomials are positive-definite in the range [0,1].
void fillBuffer() const
RooAbsRealLValue const & x() const
RooArgList const & coefList() const
double xmax() const
double xmin() const
Bifurcated Gaussian p.d.f with different widths on left and right side of maximum value.
RooAbsReal const & getSigmaL() const
Get the left sigma parameter.
RooAbsReal const & getSigmaR() const
Get the right sigma parameter.
RooAbsReal const & getX() const
Get the x variable.
RooAbsReal const & getMean() const
Get the mean parameter.
PDF implementing the Crystal Ball line shape.
Definition RooCBShape.h:24
RooAbsReal const & getSigma() const
Definition RooCBShape.h:43
RooAbsReal const & getM() const
Definition RooCBShape.h:41
RooAbsReal const & getN() const
Definition RooCBShape.h:45
RooAbsReal const & getM0() const
Definition RooCBShape.h:42
RooAbsReal const & getAlpha() const
Definition RooCBShape.h:44
Object to represent discrete states.
Definition RooCategory.h:28
value_type getCurrentIndex() const final
Return current index.
Definition RooCategory.h:40
Chebychev polynomial p.d.f.
RooAbsReal const & x() const
RooArgList const & coefList() const
const char * refRangeName() const
Represents a constant real-valued object.
Definition RooConstVar.h:23
Calculates the sum of the -(log) likelihoods of a set of RooAbsPfs that represent constraint function...
const RooArgList & list()
Container class to hold N-dimensional binned data.
Definition RooDataHist.h:40
double sum(bool correctForBinSize, bool inverseCorr=false) const
Return the sum of the weights of all bins in the histogram.
std::vector< std::unique_ptr< const RooAbsBinning > > const & getBinnings() const
std::string declWeightArrayForCodeSquash(RooFit::Experimental::CodegenContext &ctx, bool correctForBinSize) const
double weight(std::size_t i) const
Return weight of i-th bin.
std::string calculateTreeIndexForCodeSquash(RooFit::Experimental::CodegenContext &ctx, const RooAbsCollection &coords, bool reverse=false) const
The class RooEffProd implements the product of a PDF with an efficiency function.
Definition RooEffProd.h:19
RooAbsReal const & pdf() const
Definition RooEffProd.h:31
RooAbsReal const & eff() const
Definition RooEffProd.h:32
A PDF helper class to fit efficiencies parameterized by a supplied function F.
RooAbsCategory const & cat() const
RooAbsReal const & effFunc() const
std::string sigCatName() const
Exponential PDF.
bool negateCoefficient() const
RooAbsReal const & coefficient() const
Get the coefficient "c".
RooAbsReal const & variable() const
Get the x variable.
RooExtendPdf is a wrapper around an existing PDF that adds a parameteric extended likelihood term to ...
RooAbsPdf const & pdf() const
A RooProdPdf with a fixed normalization set can be replaced by this class.
Definition RooProdPdf.h:212
RooArgSet const * partList() const
Definition RooProdPdf.h:267
RooAbsReal const * rearrangedDen() const
Definition RooProdPdf.h:262
RooAbsReal const * rearrangedNum() const
Definition RooProdPdf.h:258
RooAbsReal const & normIntegral() const
RooAbsPdf const & pdf() const
A class to maintain the context for squashing of RooFit models into code.
void addToGlobalScope(std::string const &str)
Adds the given string to the string block that will be emitted at the top of the squashed function.
std::string const & getResult(RooAbsArg const &arg)
Gets the result for the given node using the node name.
std::string getTmpVarName() const
Get a unique variable name to be used in the generated code.
void addResult(RooAbsArg const *key, std::string const &value)
A function to save an expression that includes/depends on the result of the input node.
void addToCodeBody(RooAbsArg const *klass, std::string const &in)
Adds the input string to the squashed code body.
std::unique_ptr< LoopScope > beginLoop(RooAbsArg const *in)
Create a RAII scope for iterating over vector observables.
void collectFunction(std::string const &name)
Register a function that is only know to the interpreter to the context.
void addVecObs(const char *key, int idx)
Since the squashed code represents all observables as a single flattened array, it is important to ke...
int observableIndexOf(const RooAbsArg &arg) const
std::string buildArg(RooAbsCollection const &x, std::string const &arrayType="double")
Function to save a RooListProxy as an array in the squashed code.
std::string buildFunction(RooAbsArg const &arg, std::unordered_set< RooFit::Detail::DataKey > const &dependsOnData={})
Assemble and return the final code with the return expression and global statements.
std::string buildCall(std::string const &funcname, Args_t const &...args)
Build the code to call the function with name funcname, passing some arguments.
A RooFormulaVar is a generic implementation of a real-valued object, which takes a RooArgList of serv...
const RooArgList & dependents() const
std::string getUniqueFuncName() const
RooCFunction1Binding is a templated implementation of class RooAbsReal that binds generic C(++) funct...
RooAbsReal const & variable() const
RooAbsReal const & variable() const
RooFunctorBinding makes math functions from ROOT usable in RooFit.
RooArgList const & variables() const
RooFunctorPdfBinding makes math functions from ROOT usable as PDFs in RooFit.
RooArgList const & variables() const
Implementation of the Gamma PDF for RooFit/RooStats.
Definition RooGamma.h:20
RooAbsReal const & getX() const
Definition RooGamma.h:34
RooAbsReal const & getGamma() const
Definition RooGamma.h:35
RooAbsReal const & getBeta() const
Definition RooGamma.h:36
RooAbsReal const & getMu() const
Definition RooGamma.h:37
Plain Gaussian p.d.f.
Definition RooGaussian.h:24
RooAbsReal const & getX() const
Get the x variable.
Definition RooGaussian.h:45
RooAbsReal const & getMean() const
Get the mean parameter.
Definition RooGaussian.h:48
RooAbsReal const & getSigma() const
Get the sigma parameter.
Definition RooGaussian.h:51
Implementation of a probability density function that takes a RooArgList of servers and a C++ express...
const RooArgList & dependents() const
std::string getUniqueFuncName() const
A real-valued function sampled from a multidimensional histogram.
Definition RooHistFunc.h:31
Int_t getInterpolationOrder() const
Return histogram interpolation order.
Definition RooHistFunc.h:69
bool getCdfBoundaries() const
If true, special boundary conditions for c.d.f.s are used.
Definition RooHistFunc.h:85
RooDataHist & dataHist()
Return RooDataHist that is represented.
Definition RooHistFunc.h:45
RooArgSet const & variables() const
A probability density function sampled from a multidimensional histogram.
Definition RooHistPdf.h:30
Int_t getInterpolationOrder() const
Definition RooHistPdf.h:59
RooDataHist & dataHist()
Definition RooHistPdf.h:42
bool haveUnitNorm() const
Definition RooHistPdf.h:82
bool getCdfBoundaries() const
Definition RooHistPdf.h:73
RooArgSet const & variables() const
Definition RooHistPdf.h:98
Landau distribution p.d.f.
Definition RooLandau.h:24
RooAbsReal const & getSigma() const
Definition RooLandau.h:42
RooAbsReal const & getMean() const
Definition RooLandau.h:41
RooAbsReal const & getX() const
Definition RooLandau.h:40
RooFit Lognormal PDF.
bool useStandardParametrization() const
RooAbsReal const & getMedian() const
Get the median parameter.
RooAbsReal const & getShapeK() const
Get the shape parameter.
RooAbsReal const & getX() const
Get the x variable.
The class RooMultiPdf allows for the creation of a RooMultiPdf object, which can switch between previ...
Definition RooMultiPdf.h:9
const RooCategoryProxy & indexCategory() const
Definition RooMultiPdf.h:25
RooAbsPdf * getPdf(int index) const
Definition RooMultiPdf.h:31
int getNumPdfs() const
Definition RooMultiPdf.h:23
const RooListProxy & getPdfList() const
Definition RooMultiPdf.h:26
Multivariate Gaussian p.d.f.
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Handle full integral here.
const RooArgList & xVec() const
const TMatrixDSym & covarianceMatrixInverse() const
const RooArgList & muVec() const
RooONNXFunction wraps an ONNX model as a RooAbsReal, allowing it to be used as a building block in li...
RooArgList const & inputTensorList(int iTensor) const
std::size_t nInputTensors() const
std::string outerWrapperName() const
A histogram function that assigns scale parameters to every bin.
const RooArgList & paramList() const
const RooArgList & xList() const
const RooDataHist & dataHist() const
bool relParam() const
Poisson pdf.
Definition RooPoisson.h:19
RooAbsReal const & getX() const
Get the x variable.
Definition RooPoisson.h:45
bool getProtectNegativeMean() const
Definition RooPoisson.h:42
bool getNoRounding() const
Definition RooPoisson.h:37
RooAbsReal const & getMean() const
Get the mean parameter.
Definition RooPoisson.h:48
A RooAbsReal implementing a polynomial in terms of a list of RooAbsReal coefficients.
Definition RooPolyVar.h:25
RooArgList const & coefList() const
Definition RooPolyVar.h:38
int lowestOrder() const
Definition RooPolyVar.h:39
RooAbsReal const & x() const
Definition RooPolyVar.h:37
RooPolynomial implements a polynomial p.d.f of the form.
RooAbsReal const & x() const
Get the x variable.
int lowestOrder() const
Return the order for the first coefficient in the list.
RooArgList const & coefList() const
Get the coefficient list.
Represents the product of a given set of RooAbsReal objects.
Definition RooProduct.h:29
const RooArgList & realComponents() const
Definition RooProduct.h:50
Represents the ratio of two RooAbsReal objects.
Definition RooRatio.h:21
RooAbsReal const & numerator() const
Definition RooRatio.h:34
RooAbsReal const & denominator() const
Definition RooRatio.h:35
Performs hybrid numerical/analytical integrals of RooAbsReal objects.
const RooArgSet & numIntRealVars() const
RooArgSet intVars() const
const RooAbsReal & integrand() const
const RooArgSet & numIntCatVars() const
const char * intRange() const
const RooArgList & coefList() const
const RooArgList & funcList() const
Implements a PDF constructed from a sum of functions:
const RooArgList & funcList() const
const RooArgList & funcIntListFromCache(Int_t code, const char *rangeName=nullptr) const
Collect the list of functions to be integrated from the cache.
const RooArgList & coefList() const
Variable that can be changed from the outside.
Definition RooRealVar.h:37
A RooAbsReal implementation that calculates the plain fraction of sum of RooAddPdf components from a ...
RooArgList const & variables() const
const std::vector< int > & interpolationCodes() const
const RooListProxy & variables() const
const std::vector< double > & high() const
const std::vector< double > & low() const
Flat p.d.f.
Definition RooUniform.h:24
double analyticalIntegral(Int_t code, const char *rangeName=nullptr) const override
Implement analytical integral.
The RooWrapperPdf is a class that can be used to convert a function into a PDF.
RooAbsReal const & function() const
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:84
const char * GetName() const override
Returns name of object.
Definition TNamed.h:49
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:224
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
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
std::string makeValidVarName(std::string const &in)
void codegenImpl(RooFit::Detail::RooFixedProdPdf &arg, CodegenContext &ctx)
void functorCodegenImpl(RooArg_t &arg, RooArgList const &variables, CodegenContext &ctx)
std::string codegenIntegralImpl(RooAbsReal &arg, int code, const char *rangeName, CodegenContext &ctx)
This function defines the analytical integral translation for the class.
@ InputArguments
TMarker m
Definition textangle.C:8
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2338