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 <RooParamHistFunc.h>
46#include <RooPoisson.h>
47#include <RooPolyVar.h>
48#include <RooPolynomial.h>
49#include <RooProdPdf.h>
50#include <RooProduct.h>
51#include <RooRatio.h>
52#include <RooRealIntegral.h>
53#include <RooRealSumFunc.h>
54#include <RooRealSumPdf.h>
55#include <RooRealVar.h>
60#include <RooUniform.h>
61#include <RooWrapperPdf.h>
62
63#include "RooFitImplHelpers.h"
64
65#include <TInterpreter.h>
66
67#include <unordered_set>
68
69namespace RooFit::Experimental {
70
71namespace {
72
73// Return a stringy-field version of the value, formatted to maximum precision.
74std::string doubleToString(double val)
75{
76 std::stringstream ss;
77 ss << std::setprecision(std::numeric_limits<double>::max_digits10) << val;
78 return ss.str();
79}
80
81std::string mathFunc(std::string const &name)
82{
83 return "RooFit::Detail::MathFuncs::" + name;
84}
85
86void rooHistTranslateImpl(RooAbsArg const &arg, CodegenContext &ctx, int intOrder, RooDataHist const &dataHist,
87 const RooArgSet &obs, bool correctForBinSize, bool cdfBoundaries)
88{
89 if (intOrder != 0 && !(!cdfBoundaries && !correctForBinSize && intOrder == 1 && obs.size() == 1)) {
90 ooccoutE(&arg, InputArguments) << "RooHistPdf::weight(" << arg.GetName()
91 << ") ERROR: codegen currently only supports non-interpolation cases."
92 << std::endl;
93 return;
94 }
95
96 if (intOrder == 1) {
97 RooAbsBinning const &binning = *dataHist.getBinnings()[0];
98 std::string weightArr = dataHist.declWeightArrayForCodeSquash(ctx, correctForBinSize);
99 ctx.addResult(&arg, ctx.buildCall(mathFunc("interpolate1d"), binning.lowBound(), binning.highBound(), *obs[0],
100 binning.numBins(), weightArr));
101 return;
102 }
103 std::string const &offset = dataHist.calculateTreeIndexForCodeSquash(ctx, obs);
104 std::string weightArr = dataHist.declWeightArrayForCodeSquash(ctx, correctForBinSize);
105 ctx.addResult(&arg, "*(" + weightArr + " + " + offset + ")");
106}
107
108std::string realSumPdfTranslateImpl(CodegenContext &ctx, RooAbsArg const &arg, RooArgList const &funcList,
109 RooArgList const &coefList, bool normalize)
110{
111 bool noLastCoeff = funcList.size() != coefList.size();
112
113 std::string const &funcName = ctx.buildArg(funcList);
114 std::string const &coeffName = ctx.buildArg(coefList);
115 std::string const &coeffSize = std::to_string(coefList.size());
116
117 std::string sum = ctx.getTmpVarName();
118 std::string coeffSum = ctx.getTmpVarName();
119 ctx.addToCodeBody(&arg, "double " + sum + " = 0;\ndouble " + coeffSum + "= 0;\n");
120
121 std::string iterator = "i_" + ctx.getTmpVarName();
122 std::string subscriptExpr = "[" + iterator + "]";
123
124 std::string code = "for(int " + iterator + " = 0; " + iterator + " < " + coeffSize + "; " + iterator + "++) {\n" +
125 sum + " += " + funcName + subscriptExpr + " * " + coeffName + subscriptExpr + ";\n";
126 code += coeffSum + " += " + coeffName + subscriptExpr + ";\n";
127 code += "}\n";
128
129 if (noLastCoeff) {
130 code += sum + " += " + funcName + "[" + coeffSize + "]" + " * (1 - " + coeffSum + ");\n";
131 } else if (normalize) {
132 code += sum + " /= " + coeffSum + ";\n";
133 }
134 ctx.addToCodeBody(&arg, code);
135
136 return sum;
137}
138
139} // namespace
140
142{
143 if (arg.isRearranged()) {
144 ctx.addResult(&arg, ctx.buildCall(mathFunc("ratio"), *arg.rearrangedNum(), *arg.rearrangedDen()));
145 } else {
146 ctx.addResult(&arg, ctx.buildCall(mathFunc("product"), *arg.partList(), arg.partList()->size()));
147 }
148}
149
151{
152 std::string const &idx = arg.dataHist().calculateTreeIndexForCodeSquash(ctx, arg.dataVars(), true);
153 std::string const &paramNames = ctx.buildArg(arg.paramList());
154
155 ctx.addResult(&arg, paramNames + "[" + idx + "]");
156}
157
159{
160 auto const &interpCodes = arg.interpolationCodes();
161
162 std::size_t n = interpCodes.size();
163
164 std::string resName = "total_" + ctx.getTmpVarName();
165 for (std::size_t i = 0; i < n; ++i) {
166 if (interpCodes[i] != interpCodes[0]) {
168 << "FlexibleInterpVar::evaluate ERROR: Code Squashing AD does not yet support having "
169 "different interpolation codes for the same class object "
170 << std::endl;
171 }
172 }
173
174 // The PiecewiseInterpolation class is used in the context of HistFactory
175 // models, where is is always used the same way: all RooAbsReals in _lowSet,
176 // _histSet, and also nominal are 1D RooHistFuncs with with same structure.
177 //
178 // Therefore, we can make a big optimization: we get the bin index only once
179 // here in the generated code for PiecewiseInterpolation. Then, we also
180 // rearrange the histogram data in such a way that we can always pass the
181 // same arrays to the free function that implements the interpolation, just
182 // with a dynamic offset calculated from the bin index.
183 RooDataHist const &nomHist = dynamic_cast<RooHistFunc const &>(*arg.nominalHist()).dataHist();
184 int nBins = nomHist.numEntries();
185 std::vector<double> valsNominal;
186 std::vector<double> valsLow;
187 std::vector<double> valsHigh;
188 for (int i = 0; i < nBins; ++i) {
189 valsNominal.push_back(nomHist.weight(i));
190 }
191 for (int i = 0; i < nBins; ++i) {
192 for (std::size_t iParam = 0; iParam < n; ++iParam) {
193 valsLow.push_back(dynamic_cast<RooHistFunc const &>(arg.lowList()[iParam]).dataHist().weight(i));
194 valsHigh.push_back(dynamic_cast<RooHistFunc const &>(arg.highList()[iParam]).dataHist().weight(i));
195 }
196 }
197 std::string idxName = ctx.getTmpVarName();
198 std::string valsNominalStr = ctx.buildArg(valsNominal);
199 std::string valsLowStr = ctx.buildArg(valsLow);
200 std::string valsHighStr = ctx.buildArg(valsHigh);
201 std::string nStr = std::to_string(n);
202 std::string code;
203
204 std::string lowName = ctx.getTmpVarName();
205 std::string highName = ctx.getTmpVarName();
206 std::string nominalName = ctx.getTmpVarName();
207 code +=
208 "unsigned int " + idxName + " = " +
209 nomHist.calculateTreeIndexForCodeSquash(ctx, dynamic_cast<RooHistFunc const &>(*arg.nominalHist()).variables()) +
210 ";\n";
211 code += "double const* " + lowName + " = " + valsLowStr + " + " + nStr + " * " + idxName + ";\n";
212 code += "double const* " + highName + " = " + valsHighStr + " + " + nStr + " * " + idxName + ";\n";
213 code += "double " + nominalName + " = *(" + valsNominalStr + " + " + idxName + ");\n";
214
215 std::string funcCall = ctx.buildCall(mathFunc("flexibleInterp"), interpCodes[0], arg.paramList(), n, lowName,
216 highName, 1.0, nominalName, 0.0);
217 code += "double " + resName + " = " + funcCall + ";\n";
218
219 if (arg.positiveDefinite()) {
220 code += resName + " = " + resName + " < 0 ? 0 : " + resName + ";\n";
221 }
222
223 ctx.addToCodeBody(&arg, code);
224 ctx.addResult(&arg, resName);
225}
226
227////////////////////////////////////////////////////////////////////////////////
228/// This function defines a translation for each RooAbsReal based object that can be used
229/// to express the class as simple C++ code. The function adds the code represented by
230/// each class as an std::string (that is later concatenated with code strings from translate calls)
231/// to form the C++ code that AD tools can understand. Any class that wants to support AD, has to
232/// implement this function.
233///
234/// \param[in] ctx An object to manage auxiliary information for code-squashing. Also takes the
235/// code string that this class outputs into the squashed code through the 'addToCodeBody' function.
237{
238 std::stringstream errorMsg;
239 errorMsg << "Translate function for class \"" << arg.ClassName() << "\" has not yet been implemented.";
240 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
241 return ctx.addResult(&arg, "1.0");
242}
243
245{
246 ctx.addResult(&arg, realSumPdfTranslateImpl(ctx, arg, arg.pdfList(), arg.coefList(), true));
247}
248
250{
251 auto const &covI = arg.covarianceMatrixInverse();
252 std::span<const double> covISpan{covI.GetMatrixArray(), static_cast<size_t>(covI.GetNoElements())};
253 ctx.addResult(&arg,
254 ctx.buildCall(mathFunc("multiVarGaussian"), arg.xVec().size(), arg.xVec(), arg.muVec(), covISpan));
255}
256
258{
259 int numPdfs = arg.getNumPdfs();
260
261 // MathFunc call
262
263 // The value of this number should be discussed. Beyound a certain number of
264 // indices MathFunc call becomes more efficient.
265 if (numPdfs > 2) {
266 ctx.addResult(&arg, ctx.buildCall(mathFunc("multipdf"), arg.indexCategory(), arg.getPdfList()));
267
268 std::cout << "MathFunc call used\n";
269
270 } else {
271
272 // Ternary nested expression
273 std::string indexExpr = ctx.getResult(arg.indexCategory());
274
275 // int numPdfs = arg.getNumPdfs();
276 std::string expr;
277
278 for (int i = 0; i < numPdfs; ++i) {
279 RooAbsPdf *pdf = arg.getPdf(i);
280 std::string pdfExpr = ctx.getResult(*pdf);
281
282 expr += "(" + indexExpr + " == " + std::to_string(i) + " ? (" + pdfExpr + ") : ";
283 }
284
285 expr += "0.0";
286 expr += std::string(numPdfs, ')'); // Close all ternary operators
287
288 ctx.addResult(&arg, expr);
289 std::cout << "Ternary expression call used \n";
290 }
291}
292
293// RooCategory index added.
295{
296 int idx = ctx.observableIndexOf(arg);
297 if (idx < 0) {
298
299 idx = 1;
300 ctx.addVecObs(arg.GetName(), idx);
301 }
302
303 std::string result = std::to_string(arg.getCurrentIndex());
304 ctx.addResult(&arg, result);
305}
306
308{
309 if (arg.list().empty()) {
310 ctx.addResult(&arg, "0.0");
311 }
312 std::string result;
313 if (arg.list().size() > 1)
314 result += "(";
315
316 std::size_t i = 0;
317 for (auto *component : static_range_cast<RooAbsReal *>(arg.list())) {
318
319 if (!dynamic_cast<RooFit::Detail::RooNLLVarNew *>(component) || arg.list().size() == 1) {
320 result += ctx.getResult(*component);
321 ++i;
322 if (i < arg.list().size())
323 result += '+';
324 continue;
325 }
326 result += ctx.buildFunction(*component, ctx.dependsOnData()) + "(params, obs, xlArr)";
327 ++i;
328 if (i < arg.list().size())
329 result += '+';
330 }
331 if (arg.list().size() > 1)
332 result += ')';
333 ctx.addResult(&arg, result);
334}
335
337{
338 arg.fillBuffer();
339 ctx.addResult(&arg, ctx.buildCall(mathFunc("bernstein"), arg.x(), arg.xmin(), arg.xmax(), arg.coefList(),
340 arg.coefList().size()));
341}
342
344{
345 ctx.addResult(&arg,
346 ctx.buildCall(mathFunc("bifurGauss"), arg.getX(), arg.getMean(), arg.getSigmaL(), arg.getSigmaR()));
347}
348
350{
351 ctx.addResult(
352 &arg, ctx.buildCall(mathFunc("cbShape"), arg.getM(), arg.getM0(), arg.getSigma(), arg.getAlpha(), arg.getN()));
353}
354
356{
357 // first bring the range of the variable _x to the normalised range [-1, 1]
358 // calculate sum_k c_k T_k(x) where x is given in the normalised range,
359 // c_0 = 1, and the higher coefficients are given in _coefList
360 double xmax = static_cast<RooAbsRealLValue const &>(arg.x()).getMax(arg.refRangeName());
361 double xmin = static_cast<RooAbsRealLValue const &>(arg.x()).getMin(arg.refRangeName());
362
363 ctx.addResult(&arg,
364 ctx.buildCall(mathFunc("chebychev"), arg.coefList(), arg.coefList().size(), arg.x(), xmin, xmax));
365}
366
368{
369 ctx.addResult(&arg, doubleToString(arg.getVal()));
370}
371
373{
374 ctx.addResult(&arg, ctx.buildCall(mathFunc("constraintSum"), arg.list(), arg.list().size()));
375}
376
377// Generate RooFit codegen wrappers for RooFunctorBinding and similar objects,
378// emitting both the primal function call and its gradient pullback for
379// Clad-based AD.
380template <class RooArg_t>
381void functorCodegenImpl(RooArg_t &arg, RooArgList const &variables, CodegenContext &ctx)
382{
383 if (!arg.function()->HasGradient()) {
384 std::stringstream errorMsg;
385 errorMsg << "Functor wrapped by \"" << arg.GetName() << "\" doesn't provide a gradient function."
386 << " RooFit codegen is therefore not supported.";
387 oocoutE(&arg, InputArguments) << errorMsg.str() << std::endl;
388 throw std::runtime_error(errorMsg.str());
389 }
390
391 std::string funcAddrStr = TString::Format("0x%zx", reinterpret_cast<std::size_t>(arg.function())).Data();
392 std::string wrapperName = "roo_functor_" + funcAddrStr;
393
394 static std::unordered_set<std::string> wrapperNames;
395
396 if (wrapperNames.find(wrapperName) == wrapperNames.end()) {
397
399
400 std::string pullbackName = wrapperName + "_pullback";
401 std::string nStr = std::to_string(std::size(variables));
402
403 std::string type;
404 if constexpr (std::is_same_v<RooArg_t, RooFunctor1DBinding> || std::is_same_v<RooArg_t, RooFunctor1DPdfBinding>)
405 type = "::ROOT::Math::IGradientFunctionOneDim";
406 else
407 type = "::ROOT::Math::IGradientFunctionMultiDim";
408
409 std::string funcAddrCasted = "reinterpret_cast<" + type + " const *>(" + funcAddrStr + ")";
410
411 std::string code;
412
413 code += "double " + wrapperName +
414 "(double const *x) {\n"
415 " return " +
417 "->operator()(x);\n"
418 "}\n\n"
419 "namespace clad::custom_derivatives {\n\n"
420 "void " +
422 "(double const* x, double d_y, double *d_x) {\n"
423 " double output[" +
424 nStr +
425 "]{};\n"
426 " " +
428 "->Gradient(x, output);\n"
429 " for (int i = 0; i < " +
430 nStr +
431 "; ++i) {\n"
432 " d_x[i] += output[i] * d_y;\n"
433 " }\n"
434 "}\n"
435 "} // namespace clad::custom_derivatives\n";
436
437 gInterpreter->Declare(code.c_str());
438 }
439
440 ctx.addResult(&arg, ctx.buildCall(wrapperName, variables));
441}
442
444{
445 functorCodegenImpl(arg, arg.variable(), ctx);
446}
447
449{
450 functorCodegenImpl(arg, arg.variable(), ctx);
451}
452
454{
455 functorCodegenImpl(arg, arg.variables(), ctx);
456}
457
459{
460 functorCodegenImpl(arg, arg.variables(), ctx);
461}
462
464{
465 ctx.addResult(&arg, ctx.buildCall("TMath::GammaDist", arg.getX(), arg.getGamma(), arg.getMu(), arg.getBeta()));
466}
467
469{
470 arg.getVal(); // to trigger the creation of the TFormula
471 std::string funcName = arg.getUniqueFuncName();
473 // We have to force the array type to be "double" because that's what the
474 // declared function wrapped by the TFormula expects.
475 auto inputVar = ctx.buildArg(arg.dependents(), /*arrayType=*/"double");
476 ctx.addResult(&arg, funcName + "(" + inputVar + ")");
477}
478
480{
481 ctx.addResult(&arg, ctx.buildCall(mathFunc("effProd"), arg.eff(), arg.pdf()));
482}
483
485{
486 RooAbsCategory const &cat = arg.cat();
487 int sigCatIndex = cat.lookupIndex(arg.sigCatName());
488 ctx.addResult(&arg, ctx.buildCall(mathFunc("efficiency"), arg.effFunc(), cat, sigCatIndex));
489}
490
492{
493 // Build a call to the stateless exponential defined later.
494 std::string coef;
495 if (arg.negateCoefficient()) {
496 coef += "-";
497 }
498 coef += ctx.getResult(arg.coefficient());
499 ctx.addResult(&arg, "std::exp(" + coef + " * " + ctx.getResult(arg.variable()) + ")");
500}
501
503{
504 // Use the result of the underlying pdf.
505 ctx.addResult(&arg, ctx.getResult(arg.pdf()));
506}
507
509{
510 // Build a call to the stateless gaussian defined later.
511 ctx.addResult(&arg, ctx.buildCall(mathFunc("gaussian"), arg.getX(), arg.getMean(), arg.getSigma()));
512}
513
515{
516 arg.getVal(); // to trigger the creation of the TFormula
517 std::string funcName = arg.getUniqueFuncName();
519 // We have to force the array type to be "double" because that's what the
520 // declared function wrapped by the TFormula expects.
521 auto inputVar = ctx.buildArg(arg.dependents(), /*arrayType=*/"double");
522 ctx.addResult(&arg, funcName + "(" + inputVar + ")");
523}
524
526{
527 rooHistTranslateImpl(arg, ctx, arg.getInterpolationOrder(), arg.dataHist(), arg.variables(), false,
528 arg.getCdfBoundaries());
529}
530
532{
533 rooHistTranslateImpl(arg, ctx, arg.getInterpolationOrder(), arg.dataHist(), arg.variables(), !arg.haveUnitNorm(),
534 arg.getCdfBoundaries());
535}
536
538{
539 ctx.addResult(&arg, ctx.buildCall(mathFunc("landau"), arg.getX(), arg.getMean(), arg.getSigma()));
540}
541
543{
544 std::string funcName = arg.useStandardParametrization() ? "logNormalEvaluateStandard" : "logNormal";
545 ctx.addResult(&arg, ctx.buildCall(mathFunc(funcName), arg.getX(), arg.getShapeK(), arg.getMedian()));
546}
547
548void codegenImpl(RooFit::Detail::RooNLLVarNew &arg, CodegenContext &ctx)
549{
550 if (arg.binnedL() && !arg.pdf().getAttribute("BinnedLikelihoodActiveYields")) {
551 std::stringstream errorMsg;
552 errorMsg << "codegen: binned likelihood optimization is only supported when raw pdf "
553 "values can be interpreted as yields."
554 << " This is not the case for HistFactory models written with ROOT versions before 6.26.00";
555 oocoutE(&arg, InputArguments) << errorMsg.str() << std::endl;
556 throw std::runtime_error(errorMsg.str());
557 }
558
559 std::string weightSumName = RooFit::Detail::makeValidVarName(arg.GetName()) + "WeightSum";
560 std::string resName = RooFit::Detail::makeValidVarName(arg.GetName()) + "Result";
561 ctx.addResult(&arg, resName);
562 ctx.addToGlobalScope("double " + weightSumName + " = 0.0;\n");
563 ctx.addToGlobalScope("double " + resName + " = 0.0;\n");
564
565 const bool needWeightSum = arg.expectedEvents() || arg.simCount() > 1;
566
567 if (needWeightSum) {
568 auto scope = ctx.beginLoop(&arg);
569 ctx.addToCodeBody(weightSumName + " += " + ctx.getResult(arg.weightVar()) + ";\n");
570 }
571 if (arg.simCount() > 1) {
572 std::string simCountStr = std::to_string(static_cast<double>(arg.simCount()));
573 ctx.addToCodeBody(resName + " += " + weightSumName + " * std::log(" + simCountStr + ");\n");
574 }
575
576 // Begin loop scope for the observables and weight variable. If the weight
577 // is a scalar, the context will ignore it for the loop scope. The closing
578 // brackets of the loop is written at the end of the scopes lifetime.
579 {
580 auto scope = ctx.beginLoop(&arg);
581 std::string term = ctx.buildCall(mathFunc("nll"), arg.pdf(), arg.weightVar(), arg.binnedL(), 0);
582 ctx.addToCodeBody(&arg, resName + " += " + term + ";");
583 }
584 if (arg.expectedEvents()) {
585 std::string expected = ctx.getResult(*arg.expectedEvents());
586 ctx.addToCodeBody(resName + " += " + expected + " - " + weightSumName + " * std::log(" + expected + ");\n");
587 }
588}
589
591{
592 // For now just return function/normalization integral.
593 ctx.addResult(&arg, ctx.getResult(arg.pdf()) + "/" + ctx.getResult(arg.normIntegral()));
594}
595
597{
598 std::string const &idx = arg.dataHist().calculateTreeIndexForCodeSquash(ctx, arg.xList());
599 std::string arrName = ctx.buildArg(arg.paramList());
600 std::stringstream result;
601 result << arrName << "[" << idx << "]";
602 if (arg.relParam()) {
603 // get weight[idx] * binv[idx]. Here we get the bin volume for the first element as we assume the distribution to
604 // be binned uniformly.
605 double binV = arg.dataHist().binVolume(0);
606 std::string weightArr = arg.dataHist().declWeightArrayForCodeSquash(ctx, false);
607 result << " * *(" << weightArr << " + " << idx + ") * " << doubleToString(binV);
608 }
609 ctx.addResult(&arg, result.str());
610}
611
613{
614 std::string xName = ctx.getResult(arg.getX());
615 if (!arg.getNoRounding())
616 xName = "std::floor(" + xName + ")";
617
618 ctx.addResult(&arg, ctx.buildCall(mathFunc("poisson"), xName, arg.getMean()));
619}
620
622{
623 const unsigned sz = arg.coefList().size();
624 if (!sz) {
625 ctx.addResult(&arg, std::to_string(arg.lowestOrder() ? 1. : 0.));
626 return;
627 }
628
629 ctx.addResult(&arg, ctx.buildCall(mathFunc("polynomial"), arg.coefList(), sz, arg.lowestOrder(), arg.x()));
630}
631
633{
634 const unsigned sz = arg.coefList().size();
635 if (!sz) {
636 ctx.addResult(&arg, std::to_string(arg.lowestOrder() ? 1. : 0.));
637 return;
638 }
639
640 ctx.addResult(&arg, ctx.buildCall(mathFunc("polynomial<true>"), arg.coefList(), sz, arg.lowestOrder(), arg.x()));
641}
642
644{
645 ctx.addResult(&arg, ctx.buildCall(mathFunc("product"), arg.realComponents(), arg.realComponents().size()));
646}
647
649{
650 ctx.addResult(&arg, ctx.buildCall(mathFunc("ratio"), arg.numerator(), arg.denominator()));
651}
652
653namespace {
654
655std::string codegenIntegral(RooAbsReal &arg, int code, const char *rangeName, CodegenContext &ctx)
656{
657 using Func = std::string (*)(RooAbsReal &, int, const char *, CodegenContext &);
658
659 Func func;
660
661 TClass *tclass = arg.IsA();
662
663 // Cache the overload resolutions
664 static std::unordered_map<TClass *, Func> dispatchMap;
665
666 auto found = dispatchMap.find(tclass);
667
668 if (found != dispatchMap.end()) {
669 func = found->second;
670 } else {
671 // Can probably done with CppInterop in the future to avoid string manipulation.
672 std::stringstream cmd;
673 cmd << "&RooFit::Experimental::CodegenIntegralImplCaller<" << tclass->GetName() << ">::call;";
674 func = reinterpret_cast<Func>(gInterpreter->ProcessLine(cmd.str().c_str()));
675 dispatchMap[tclass] = func;
676 }
677
678 return func(arg, code, rangeName, ctx);
679}
680
681} // namespace
682
684{
685 if (arg.numIntCatVars().empty() && arg.numIntRealVars().empty()) {
686 ctx.addResult(&arg, codegenIntegral(const_cast<RooAbsReal &>(arg.integrand()), arg.mode(), arg.intRange(), ctx));
687 return;
688 }
689
690 if (arg.intVars().size() != 1 || arg.numIntRealVars().size() != 1) {
691 std::stringstream errorMsg;
692 errorMsg << "Only analytical integrals and 1D numeric integrals are supported for AD for class"
693 << arg.integrand().GetName();
694 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
695 throw std::runtime_error(errorMsg.str().c_str());
696 }
697
698 auto &intVar = static_cast<RooAbsRealLValue &>(*arg.numIntRealVars()[0]);
699
700 std::string obsName = ctx.getTmpVarName();
701 std::string oldIntVarResult = ctx.getResult(intVar);
702 ctx.addResult(&intVar, "obs[0]");
703
704 std::string funcName = ctx.buildFunction(arg.integrand(), {});
705
706 std::stringstream ss;
707
708 ss << "double " << obsName << "[1];\n";
709
710 std::string resName = RooFit::Detail::makeValidVarName(arg.GetName()) + "Result";
711 ctx.addResult(&arg, resName);
712 ctx.addToGlobalScope("double " + resName + " = 0.0;\n");
713
714 // TODO: once Clad has support for higher-order functions (follow also the
715 // Clad issue #637), we could refactor this code into an actual function
716 // instead of hardcoding it here as a string.
717 ss << "{\n"
718 << " const int n = 1000; // number of sampling points\n"
719 << " double d = " << intVar.getMax(arg.intRange()) << " - " << intVar.getMin(arg.intRange()) << ";\n"
720 << " double eps = d / n;\n"
721 << " #pragma clad checkpoint loop\n"
722 << " for (int i = 0; i < n; ++i) {\n"
723 << " " << obsName << "[0] = " << intVar.getMin(arg.intRange()) << " + eps * i;\n"
724 << " double tmpA = " << funcName << "(params, " << obsName << ", xlArr);\n"
725 << " " << obsName << "[0] = " << intVar.getMin(arg.intRange()) << " + eps * (i + 1);\n"
726 << " double tmpB = " << funcName << "(params, " << obsName << ", xlArr);\n"
727 << " " << resName << " += (tmpA + tmpB) * 0.5 * eps;\n"
728 << " }\n"
729 << "}\n";
730
731 ctx.addToGlobalScope(ss.str());
732
734}
735
737{
738 ctx.addResult(&arg, realSumPdfTranslateImpl(ctx, arg, arg.funcList(), arg.coefList(), false));
739}
740
742{
743 ctx.addResult(&arg, realSumPdfTranslateImpl(ctx, arg, arg.funcList(), arg.coefList(), false));
744}
745
747{
748 if (!arg.isConstant()) {
749 ctx.addResult(&arg, arg.GetName());
750 }
751 ctx.addResult(&arg, doubleToString(arg.getVal()));
752}
753
755{
756 ctx.addResult(&arg, ctx.buildCall(mathFunc("recursiveFraction"), arg.variables(), arg.variables().size()));
757}
758
760{
761 auto const &interpCodes = arg.interpolationCodes();
762
763 unsigned int n = interpCodes.size();
764
765 int interpCode = interpCodes[0];
766 // To get consistent codes with the PiecewiseInterpolation
767 if (interpCode == 4) {
768 interpCode = 5;
769 }
770
771 for (unsigned int i = 1; i < n; i++) {
772 if (interpCodes[i] != interpCodes[0]) {
774 << "FlexibleInterpVar::evaluate ERROR: Code Squashing AD does not yet support having "
775 "different interpolation codes for the same class object "
776 << std::endl;
777 }
778 }
779
780 std::string const &resName = ctx.buildCall(mathFunc("flexibleInterp"), interpCode, arg.variables(), n, arg.low(),
781 arg.high(), arg.globalBoundary(), arg.nominal(), 1.0);
782 ctx.addResult(&arg, resName);
783}
784
786{
787 ctx.addResult(&arg, "1.0");
788}
789
791{
792 ctx.addResult(&arg, ctx.getResult(arg.function()));
793}
794
795////////////////////////////////////////////////////////////////////////////////
796/// This function defines the analytical integral translation for the class.
797///
798/// \param[in] code The code that decides the integrands.
799/// \param[in] rangeName Name of the normalization range.
800/// \param[in] ctx An object to manage auxiliary information for code-squashing.
801///
802/// \returns The representative code string of the integral for the given object.
803std::string codegenIntegralImpl(RooAbsReal &arg, int, const char *, CodegenContext &)
804{
805 std::stringstream errorMsg;
806 errorMsg << "An analytical integral function for class \"" << arg.ClassName() << "\" has not yet been implemented.";
807 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
808 throw std::runtime_error(errorMsg.str().c_str());
809}
810
811std::string codegenIntegralImpl(RooBernstein &arg, int, const char *rangeName, CodegenContext &ctx)
812{
813 arg.fillBuffer(); // to get the right xmin() and xmax()
814 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
815 return ctx.buildCall(mathFunc("bernsteinIntegral"), x.getMin(rangeName), x.getMax(rangeName), arg.xmin(), arg.xmax(),
816 arg.coefList(), arg.coefList().size());
817}
818
819std::string codegenIntegralImpl(RooBifurGauss &arg, int code, const char *rangeName, CodegenContext &ctx)
820{
821 auto &constant = code == 1 ? arg.getMean() : arg.getX();
822 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(code == 1 ? arg.getX() : arg.getMean());
823
824 return ctx.buildCall(mathFunc("bifurGaussIntegral"), integrand.getMin(rangeName), integrand.getMax(rangeName),
825 constant, arg.getSigmaL(), arg.getSigmaR());
826}
827
828std::string codegenIntegralImpl(RooCBShape &arg, int /*code*/, const char *rangeName, CodegenContext &ctx)
829{
830 auto &m = dynamic_cast<RooAbsRealLValue const &>(arg.getM());
831 return ctx.buildCall(mathFunc("cbShapeIntegral"), m.getMin(rangeName), m.getMax(rangeName), arg.getM0(),
832 arg.getSigma(), arg.getAlpha(), arg.getN());
833}
834
835std::string codegenIntegralImpl(RooChebychev &arg, int, const char *rangeName, CodegenContext &ctx)
836{
837 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
838 double xmax = x.getMax(arg.refRangeName());
839 double xmin = x.getMin(arg.refRangeName());
840 unsigned int sz = arg.coefList().size();
841
842 return ctx.buildCall(mathFunc("chebychevIntegral"), arg.coefList(), sz, xmin, xmax, x.getMin(rangeName),
843 x.getMax(rangeName));
844}
845
846std::string codegenIntegralImpl(RooEfficiency &, int, const char *, CodegenContext &)
847{
848 return "1.0";
849}
850
851std::string codegenIntegralImpl(RooExponential &arg, int code, const char *rangeName, CodegenContext &ctx)
852{
853 bool isOverX = code == 1;
854
855 std::string constant;
856 if (arg.negateCoefficient() && isOverX) {
857 constant += "-";
858 }
859 constant += ctx.getResult(isOverX ? arg.coefficient() : arg.variable());
860
861 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(isOverX ? arg.variable() : arg.coefficient());
862
863 double min = integrand.getMin(rangeName);
864 double max = integrand.getMax(rangeName);
865
866 if (!isOverX && arg.negateCoefficient()) {
867 std::swap(min, max);
868 min = -min;
869 max = -max;
870 }
871
872 return ctx.buildCall(mathFunc("exponentialIntegral"), min, max, constant);
873}
874
875std::string codegenIntegralImpl(RooGamma &arg, int, const char *rangeName, CodegenContext &ctx)
876{
877 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.getX());
878 const std::string a =
879 ctx.buildCall("ROOT::Math::gamma_cdf", x.getMax(rangeName), arg.getGamma(), arg.getBeta(), arg.getMu());
880 const std::string b =
881 ctx.buildCall("ROOT::Math::gamma_cdf", x.getMin(rangeName), arg.getGamma(), arg.getBeta(), arg.getMu());
882 return a + " - " + b;
883}
884
885std::string codegenIntegralImpl(RooGaussian &arg, int code, const char *rangeName, CodegenContext &ctx)
886{
887 auto &constant = code == 1 ? arg.getMean() : arg.getX();
888 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(code == 1 ? arg.getX() : arg.getMean());
889
890 return ctx.buildCall(mathFunc("gaussianIntegral"), integrand.getMin(rangeName), integrand.getMax(rangeName),
891 constant, arg.getSigma());
892}
893
894namespace {
895
896std::string rooHistIntegralTranslateImpl(int code, RooAbsArg const &arg, RooDataHist const &dataHist,
897 const RooArgSet &obs, bool histFuncMode)
898{
899 if (((2 << obs.size()) - 1) != code) {
900 oocoutE(&arg, InputArguments) << "RooHistPdf::integral(" << arg.GetName()
901 << ") ERROR: AD currently only supports integrating over all histogram observables."
902 << std::endl;
903 return "";
904 }
905 return doubleToString(dataHist.sum(histFuncMode));
906}
907
908} // namespace
909
910std::string codegenIntegralImpl(RooHistFunc &arg, int code, const char *, CodegenContext &)
911{
912 return rooHistIntegralTranslateImpl(code, arg, arg.dataHist(), arg.variables(), true);
913}
914
915std::string codegenIntegralImpl(RooHistPdf &arg, int code, const char *, CodegenContext &)
916{
917 return rooHistIntegralTranslateImpl(code, arg, arg.dataHist(), arg.variables(), false);
918}
919
920std::string codegenIntegralImpl(RooLandau &arg, int, const char *rangeName, CodegenContext &ctx)
921{
922 // Don't do anything with "code". It can only be "1" anyway (see
923 // implementation of getAnalyticalIntegral).
924 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.getX());
925 const std::string a = ctx.buildCall("ROOT::Math::landau_cdf", x.getMax(rangeName), arg.getSigma(), arg.getMean());
926 const std::string b = ctx.buildCall("ROOT::Math::landau_cdf", x.getMin(rangeName), arg.getSigma(), arg.getMean());
927 return ctx.getResult(arg.getSigma()) + " * " + "(" + a + " - " + b + ")";
928}
929
930std::string codegenIntegralImpl(RooLognormal &arg, int, const char *rangeName, CodegenContext &ctx)
931{
932 std::string funcName = arg.useStandardParametrization() ? "logNormalIntegralStandard" : "logNormalIntegral";
933 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.getX());
934 return ctx.buildCall(mathFunc(funcName), x.getMin(rangeName), x.getMax(rangeName), arg.getMedian(), arg.getShapeK());
935}
936
937std::string codegenIntegralImpl(RooMultiVarGaussian &arg, int code, const char *rangeName, CodegenContext &)
938{
939 if (code != -1) {
940 std::stringstream errorMsg;
941 errorMsg << "Partial integrals over RooMultiVarGaussian are not supported.";
942 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
943 throw std::runtime_error(errorMsg.str().c_str());
944 }
945
947}
948
949std::string codegenIntegralImpl(RooPoisson &arg, int code, const char *rangeName, CodegenContext &ctx)
950{
951 assert(code == 1 || code == 2);
952 std::string xName = ctx.getResult(arg.getX());
953 if (!arg.getNoRounding())
954 xName = "std::floor(" + xName + ")";
955
956 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(code == 1 ? arg.getX() : arg.getMean());
957 // Since the integral function is the same for both codes, we need to make sure the indexed observables do not appear
958 // in the function if they are not required.
959 xName = code == 1 ? "0" : xName;
960 return ctx.buildCall(mathFunc("poissonIntegral"), code, arg.getMean(), xName, integrand.getMin(rangeName),
961 integrand.getMax(rangeName), arg.getProtectNegativeMean());
962}
963
964std::string codegenIntegralImpl(RooPolyVar &arg, int, const char *rangeName, CodegenContext &ctx)
965{
966 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
967 const double xmin = x.getMin(rangeName);
968 const double xmax = x.getMax(rangeName);
969 const unsigned sz = arg.coefList().size();
970 if (!sz)
971 return std::to_string(arg.lowestOrder() ? xmax - xmin : 0.0);
972
973 return ctx.buildCall(mathFunc("polynomialIntegral"), arg.coefList(), sz, arg.lowestOrder(), xmin, xmax);
974}
975
976std::string codegenIntegralImpl(RooPolynomial &arg, int, const char *rangeName, CodegenContext &ctx)
977{
978 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
979 const double xmin = x.getMin(rangeName);
980 const double xmax = x.getMax(rangeName);
981 const unsigned sz = arg.coefList().size();
982 if (!sz)
983 return std::to_string(arg.lowestOrder() ? xmax - xmin : 0.0);
984
985 return ctx.buildCall(mathFunc("polynomialIntegral<true>"), arg.coefList(), sz, arg.lowestOrder(), xmin, xmax);
986}
987
988std::string codegenIntegralImpl(RooRealSumPdf &arg, int code, const char *rangeName, CodegenContext &ctx)
989{
990 // Re-use translate, since integration is linear.
991 return realSumPdfTranslateImpl(ctx, arg, arg.funcIntListFromCache(code, rangeName), arg.coefList(), false);
992}
993
994std::string codegenIntegralImpl(RooUniform &arg, int code, const char *rangeName, CodegenContext &)
995{
996 // The integral of a uniform distribution is static, so we can just hardcode
997 // the result in a string.
999}
1000
1001} // 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:110
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
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:226
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