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 <RooChebychev.h>
24#include <RooConstVar.h>
25#include <RooConstraintSum.h>
26#include <RooEffProd.h>
27#include <RooEfficiency.h>
28#include <RooExponential.h>
29#include <RooExtendPdf.h>
32#include <RooFormulaVar.h>
33#include <RooGamma.h>
34#include <RooGaussian.h>
35#include <RooGenericPdf.h>
36#include <RooHistFunc.h>
37#include <RooHistPdf.h>
38#include <RooLandau.h>
39#include <RooLognormal.h>
40#include <RooMultiVarGaussian.h>
41#include <RooParamHistFunc.h>
42#include <RooPoisson.h>
43#include <RooPolyVar.h>
44#include <RooPolynomial.h>
45#include <RooProdPdf.h>
46#include <RooProduct.h>
47#include <RooRatio.h>
48#include <RooRealIntegral.h>
49#include <RooRealSumFunc.h>
50#include <RooRealSumPdf.h>
51#include <RooRealVar.h>
56#include <RooUniform.h>
57
58#include "RooFitImplHelpers.h"
59
60#include <TInterpreter.h>
61
62namespace RooFit {
63namespace Experimental {
64
65namespace {
66
67std::string mathFunc(std::string const &name)
68{
69 return "RooFit::Detail::MathFuncs::" + name;
70}
71
72void rooHistTranslateImpl(RooAbsArg const &arg, CodegenContext &ctx, int intOrder, RooDataHist const &dataHist,
73 const RooArgSet &obs, bool correctForBinSize, bool cdfBoundaries)
74{
75 if (intOrder != 0 && !(!cdfBoundaries && !correctForBinSize && intOrder == 1 && obs.size() == 1)) {
76 ooccoutE(&arg, InputArguments) << "RooHistPdf::weight(" << arg.GetName()
77 << ") ERROR: codegen currently only supports non-interpolation cases."
78 << std::endl;
79 return;
80 }
81
82 if (intOrder == 1) {
83 RooAbsBinning const &binning = *dataHist.getBinnings()[0];
84 std::string weightArr = dataHist.declWeightArrayForCodeSquash(ctx, correctForBinSize);
85 ctx.addResult(&arg, ctx.buildCall(mathFunc("interpolate1d"), binning.lowBound(), binning.highBound(), *obs[0],
86 binning.numBins(), weightArr));
87 return;
88 }
89 std::string const &offset = dataHist.calculateTreeIndexForCodeSquash(&arg, ctx, obs);
90 std::string weightArr = dataHist.declWeightArrayForCodeSquash(ctx, correctForBinSize);
91 ctx.addResult(&arg, "*(" + weightArr + " + " + offset + ")");
92}
93
94std::string realSumPdfTranslateImpl(CodegenContext &ctx, RooAbsArg const &arg, RooArgList const &funcList,
95 RooArgList const &coefList, bool normalize)
96{
97 bool noLastCoeff = funcList.size() != coefList.size();
98
99 std::string const &funcName = ctx.buildArg(funcList);
100 std::string const &coeffName = ctx.buildArg(coefList);
101 std::string const &coeffSize = std::to_string(coefList.size());
102
103 std::string sum = ctx.getTmpVarName();
104 std::string coeffSum = ctx.getTmpVarName();
105 ctx.addToCodeBody(&arg, "double " + sum + " = 0;\ndouble " + coeffSum + "= 0;\n");
106
107 std::string iterator = "i_" + ctx.getTmpVarName();
108 std::string subscriptExpr = "[" + iterator + "]";
109
110 std::string code = "for(int " + iterator + " = 0; " + iterator + " < " + coeffSize + "; " + iterator + "++) {\n" +
111 sum + " += " + funcName + subscriptExpr + " * " + coeffName + subscriptExpr + ";\n";
112 code += coeffSum + " += " + coeffName + subscriptExpr + ";\n";
113 code += "}\n";
114
115 if (noLastCoeff) {
116 code += sum + " += " + funcName + "[" + coeffSize + "]" + " * (1 - " + coeffSum + ");\n";
117 } else if (normalize) {
118 code += sum + " /= " + coeffSum + ";\n";
119 }
120 ctx.addToCodeBody(&arg, code);
121
122 return sum;
123}
124
125} // namespace
126
128{
129 if (arg.cache()._isRearranged) {
130 ctx.addResult(&arg, ctx.buildCall(mathFunc("ratio"), *arg.cache()._rearrangedNum, *arg.cache()._rearrangedDen));
131 } else {
132 ctx.addResult(&arg, ctx.buildCall(mathFunc("product"), arg.cache()._partList, arg.cache()._partList.size()));
133 }
134}
135
137{
138 std::string const &idx = arg.dataHist().calculateTreeIndexForCodeSquash(&arg, ctx, arg.dataVars(), true);
139 std::string const &paramNames = ctx.buildArg(arg.paramList());
140
141 ctx.addResult(&arg, paramNames + "[" + idx + "]");
142}
143
145{
146 auto const &interpCodes = arg.interpolationCodes();
147
148 std::size_t n = interpCodes.size();
149
150 std::string resName = "total_" + ctx.getTmpVarName();
151 for (std::size_t i = 0; i < n; ++i) {
152 if (interpCodes[i] != interpCodes[0]) {
154 << "FlexibleInterpVar::evaluate ERROR: Code Squashing AD does not yet support having "
155 "different interpolation codes for the same class object "
156 << std::endl;
157 }
158 }
159
160 // The PiecewiseInterpolation class is used in the context of HistFactory
161 // models, where is is always used the same way: all RooAbsReals in _lowSet,
162 // _histSet, and also nominal are 1D RooHistFuncs with with same structure.
163 //
164 // Therefore, we can make a big optimization: we get the bin index only once
165 // here in the generated code for PiecewiseInterpolation. Then, we also
166 // rearrange the histogram data in such a way that we can always pass the
167 // same arrays to the free function that implements the interpolation, just
168 // with a dynamic offset calculated from the bin index.
169 RooDataHist const &nomHist = dynamic_cast<RooHistFunc const &>(*arg.nominalHist()).dataHist();
170 int nBins = nomHist.numEntries();
171 std::vector<double> valsNominal;
172 std::vector<double> valsLow;
173 std::vector<double> valsHigh;
174 for (int i = 0; i < nBins; ++i) {
175 valsNominal.push_back(nomHist.weight(i));
176 }
177 for (int i = 0; i < nBins; ++i) {
178 for (std::size_t iParam = 0; iParam < n; ++iParam) {
179 valsLow.push_back(dynamic_cast<RooHistFunc const &>(arg.lowList()[iParam]).dataHist().weight(i));
180 valsHigh.push_back(dynamic_cast<RooHistFunc const &>(arg.highList()[iParam]).dataHist().weight(i));
181 }
182 }
183 std::string idxName = ctx.getTmpVarName();
184 std::string valsNominalStr = ctx.buildArg(valsNominal);
185 std::string valsLowStr = ctx.buildArg(valsLow);
186 std::string valsHighStr = ctx.buildArg(valsHigh);
187 std::string nStr = std::to_string(n);
188 std::string code;
189
190 std::string lowName = ctx.getTmpVarName();
191 std::string highName = ctx.getTmpVarName();
192 std::string nominalName = ctx.getTmpVarName();
193 code += "unsigned int " + idxName + " = " +
194 nomHist.calculateTreeIndexForCodeSquash(&arg, ctx,
195 dynamic_cast<RooHistFunc const &>(*arg.nominalHist()).variables()) +
196 ";\n";
197 code += "double const* " + lowName + " = " + valsLowStr + " + " + nStr + " * " + idxName + ";\n";
198 code += "double const* " + highName + " = " + valsHighStr + " + " + nStr + " * " + idxName + ";\n";
199 code += "double " + nominalName + " = *(" + valsNominalStr + " + " + idxName + ");\n";
200
201 std::string funcCall = ctx.buildCall(mathFunc("flexibleInterp"), interpCodes[0], arg.paramList(), n, lowName,
202 highName, 1.0, nominalName, 0.0);
203 code += "double " + resName + " = " + funcCall + ";\n";
204
205 if (arg.positiveDefinite()) {
206 code += resName + " = " + resName + " < 0 ? 0 : " + resName + ";\n";
207 }
208
209 ctx.addToCodeBody(&arg, code);
210 ctx.addResult(&arg, resName);
211}
212
213////////////////////////////////////////////////////////////////////////////////
214/// This function defines a translation for each RooAbsReal based object that can be used
215/// to express the class as simple C++ code. The function adds the code represented by
216/// each class as an std::string (that is later concatenated with code strings from translate calls)
217/// to form the C++ code that AD tools can understand. Any class that wants to support AD, has to
218/// implement this function.
219///
220/// \param[in] ctx An object to manage auxiliary information for code-squashing. Also takes the
221/// code string that this class outputs into the squashed code through the 'addToCodeBody' function.
223{
224 std::stringstream errorMsg;
225 errorMsg << "Translate function for class \"" << arg.ClassName() << "\" has not yet been implemented.";
226 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
227 return ctx.addResult(&arg, "1.0");
228}
229
231{
232 ctx.addResult(&arg, realSumPdfTranslateImpl(ctx, arg, arg.pdfList(), arg.coefList(), true));
233}
234
236{
237 auto const &covI = arg.covarianceMatrixInverse();
238 std::span<const double> covISpan{covI.GetMatrixArray(), static_cast<size_t>(covI.GetNoElements())};
239 ctx.addResult(&arg,
240 ctx.buildCall(mathFunc("multiVarGaussian"), arg.xVec().size(), arg.xVec(), arg.muVec(), covISpan));
241}
242
244{
245 if (arg.list().empty()) {
246 ctx.addResult(&arg, "0.0");
247 }
248 std::string result;
249 if (arg.list().size() > 1)
250 result += "(";
251
252 std::size_t i = 0;
253 for (auto *component : static_range_cast<RooAbsReal *>(arg.list())) {
254
255 if (!dynamic_cast<RooFit::Detail::RooNLLVarNew *>(component) || arg.list().size() == 1) {
256 result += ctx.getResult(*component);
257 ++i;
258 if (i < arg.list().size())
259 result += '+';
260 continue;
261 }
262 result += ctx.buildFunction(*component, ctx.outputSizes()) + "(params, obs, xlArr)";
263 ++i;
264 if (i < arg.list().size())
265 result += '+';
266 }
267 if (arg.list().size() > 1)
268 result += ')';
269 ctx.addResult(&arg, result);
270}
271
273{
274 arg.fillBuffer();
275 ctx.addResult(&arg, ctx.buildCall(mathFunc("bernstein"), arg.x(), arg.xmin(), arg.xmax(), arg.coefList(),
276 arg.coefList().size()));
277}
278
280{
281 ctx.addResult(&arg,
282 ctx.buildCall(mathFunc("bifurGauss"), arg.getX(), arg.getMean(), arg.getSigmaL(), arg.getSigmaR()));
283}
284
286{
287 ctx.addResult(
288 &arg, ctx.buildCall(mathFunc("cbShape"), arg.getM(), arg.getM0(), arg.getSigma(), arg.getAlpha(), arg.getN()));
289}
290
292{
293 // first bring the range of the variable _x to the normalised range [-1, 1]
294 // calculate sum_k c_k T_k(x) where x is given in the normalised range,
295 // c_0 = 1, and the higher coefficients are given in _coefList
296 double xmax = static_cast<RooAbsRealLValue const &>(arg.x()).getMax(arg.refRangeName());
297 double xmin = static_cast<RooAbsRealLValue const &>(arg.x()).getMin(arg.refRangeName());
298
299 ctx.addResult(&arg,
300 ctx.buildCall(mathFunc("chebychev"), arg.coefList(), arg.coefList().size(), arg.x(), xmin, xmax));
301}
302
304{
305 // Just return a stringy-field version of the const value.
306 // Formats to the maximum precision.
307 constexpr auto max_precision{std::numeric_limits<double>::digits10 + 1};
308 std::stringstream ss;
309 ss.precision(max_precision);
310 // Just use toString to make sure we do not output 'inf'.
311 // This is really ugly for large numbers...
312 ss << std::fixed << RooNumber::toString(arg.getVal());
313 ctx.addResult(&arg, ss.str());
314}
315
317{
318 ctx.addResult(&arg, ctx.buildCall(mathFunc("constraintSum"), arg.list(), arg.list().size()));
319}
320
322{
323 ctx.addResult(&arg, ctx.buildCall("TMath::GammaDist", arg.getX(), arg.getGamma(), arg.getMu(), arg.getBeta()));
324}
325
327{
328 arg.getVal(); // to trigger the creation of the TFormula
329 std::string funcName = arg.getUniqueFuncName();
330 ctx.collectFunction(funcName);
331 ctx.addResult(&arg, ctx.buildCall(funcName, arg.dependents()));
332}
333
335{
336 ctx.addResult(&arg, ctx.buildCall(mathFunc("effProd"), arg.eff(), arg.pdf()));
337}
338
340{
341 RooAbsCategory const &cat = arg.cat();
342 int sigCatIndex = cat.lookupIndex(arg.sigCatName());
343 ctx.addResult(&arg, ctx.buildCall(mathFunc("efficiency"), arg.effFunc(), cat, sigCatIndex));
344}
345
347{
348 // Build a call to the stateless exponential defined later.
349 std::string coef;
350 if (arg.negateCoefficient()) {
351 coef += "-";
352 }
353 coef += ctx.getResult(arg.coefficient());
354 ctx.addResult(&arg, "std::exp(" + coef + " * " + ctx.getResult(arg.variable()) + ")");
355}
356
358{
359 // Use the result of the underlying pdf.
360 ctx.addResult(&arg, ctx.getResult(arg.pdf()));
361}
362
364{
365 // Build a call to the stateless gaussian defined later.
366 ctx.addResult(&arg, ctx.buildCall(mathFunc("gaussian"), arg.getX(), arg.getMean(), arg.getSigma()));
367}
368
370{
371 arg.getVal(); // to trigger the creation of the TFormula
372 std::string funcName = arg.getUniqueFuncName();
373 ctx.collectFunction(funcName);
374 ctx.addResult(&arg, ctx.buildCall(funcName, arg.dependents()));
375}
376
378{
379 rooHistTranslateImpl(arg, ctx, arg.getInterpolationOrder(), arg.dataHist(), arg.variables(), false,
380 arg.getCdfBoundaries());
381}
382
384{
385 rooHistTranslateImpl(arg, ctx, arg.getInterpolationOrder(), arg.dataHist(), arg.variables(), !arg.haveUnitNorm(),
386 arg.getCdfBoundaries());
387}
388
390{
391 ctx.addResult(&arg, ctx.buildCall(mathFunc("landau"), arg.getX(), arg.getMean(), arg.getSigma()));
392}
393
395{
396 std::string funcName = arg.useStandardParametrization() ? "logNormalEvaluateStandard" : "logNormal";
397 ctx.addResult(&arg, ctx.buildCall(mathFunc(funcName), arg.getX(), arg.getShapeK(), arg.getMedian()));
398}
399
400void codegenImpl(RooFit::Detail::RooNLLVarNew &arg, CodegenContext &ctx)
401{
402 if (arg.binnedL() && !arg.pdf().getAttribute("BinnedLikelihoodActiveYields")) {
403 std::stringstream errorMsg;
404 errorMsg << "codegen: binned likelihood optimization is only supported when raw pdf "
405 "values can be interpreted as yields."
406 << " This is not the case for HistFactory models written with ROOT versions before 6.26.00";
407 oocoutE(&arg, InputArguments) << errorMsg.str() << std::endl;
408 throw std::runtime_error(errorMsg.str());
409 }
410
411 std::string weightSumName = RooFit::Detail::makeValidVarName(arg.GetName()) + "WeightSum";
412 std::string resName = RooFit::Detail::makeValidVarName(arg.GetName()) + "Result";
413 ctx.addResult(&arg, resName);
414 ctx.addToGlobalScope("double " + weightSumName + " = 0.0;\n");
415 ctx.addToGlobalScope("double " + resName + " = 0.0;\n");
416
417 const bool needWeightSum = arg.expectedEvents() || arg.simCount() > 1;
418
419 if (needWeightSum) {
420 auto scope = ctx.beginLoop(&arg);
421 ctx.addToCodeBody(weightSumName + " += " + ctx.getResult(arg.weightVar()) + ";\n");
422 }
423 if (arg.simCount() > 1) {
424 std::string simCountStr = std::to_string(static_cast<double>(arg.simCount()));
425 ctx.addToCodeBody(resName + " += " + weightSumName + " * std::log(" + simCountStr + ");\n");
426 }
427
428 // Begin loop scope for the observables and weight variable. If the weight
429 // is a scalar, the context will ignore it for the loop scope. The closing
430 // brackets of the loop is written at the end of the scopes lifetime.
431 {
432 auto scope = ctx.beginLoop(&arg);
433 std::string term = ctx.buildCall(mathFunc("nll"), arg.pdf(), arg.weightVar(), arg.binnedL(), 0);
434 ctx.addToCodeBody(&arg, resName + " += " + term + ";");
435 }
436 if (arg.expectedEvents()) {
437 std::string expected = ctx.getResult(*arg.expectedEvents());
438 ctx.addToCodeBody(resName + " += " + expected + " - " + weightSumName + " * std::log(" + expected + ");\n");
439 }
440}
441
443{
444 // For now just return function/normalization integral.
445 ctx.addResult(&arg, ctx.getResult(arg.pdf()) + "/" + ctx.getResult(arg.normIntegral()));
446}
447
449{
450 std::string const &idx = arg.dataHist().calculateTreeIndexForCodeSquash(&arg, ctx, arg.xList());
451 std::string arrName = ctx.buildArg(arg.paramList());
452 std::string result = arrName + "[" + idx + "]";
453 if (arg.relParam()) {
454 // get weight[idx] * binv[idx]. Here we get the bin volume for the first element as we assume the distribution to
455 // be binned uniformly.
456 double binV = arg.dataHist().binVolume(0);
457 std::string weightArr = arg.dataHist().declWeightArrayForCodeSquash(ctx, false);
458 result += " * *(" + weightArr + " + " + idx + ") * " + std::to_string(binV);
459 }
460 ctx.addResult(&arg, result);
461}
462
464{
465 std::string xName = ctx.getResult(arg.getX());
466 if (!arg.getNoRounding())
467 xName = "std::floor(" + xName + ")";
468
469 ctx.addResult(&arg, ctx.buildCall(mathFunc("poisson"), xName, arg.getMean()));
470}
471
473{
474 const unsigned sz = arg.coefList().size();
475 if (!sz) {
476 ctx.addResult(&arg, std::to_string(arg.lowestOrder() ? 1. : 0.));
477 return;
478 }
479
480 ctx.addResult(&arg, ctx.buildCall(mathFunc("polynomial"), arg.coefList(), sz, arg.lowestOrder(), arg.x()));
481}
482
484{
485 const unsigned sz = arg.coefList().size();
486 if (!sz) {
487 ctx.addResult(&arg, std::to_string(arg.lowestOrder() ? 1. : 0.));
488 return;
489 }
490
491 ctx.addResult(&arg, ctx.buildCall(mathFunc("polynomial<true>"), arg.coefList(), sz, arg.lowestOrder(), arg.x()));
492}
493
495{
496 ctx.addResult(&arg, ctx.buildCall(mathFunc("product"), arg.realComponents(), arg.realComponents().size()));
497}
498
500{
501 ctx.addResult(&arg, ctx.buildCall(mathFunc("ratio"), arg.numerator(), arg.denominator()));
502}
503
504namespace {
505
506std::string codegenIntegral(RooAbsReal &arg, int code, const char *rangeName, CodegenContext &ctx)
507{
508 using Func = std::string (*)(RooAbsReal &, int, const char *, CodegenContext &);
509
510 Func func;
511
512 TClass *tclass = arg.IsA();
513
514 // Cache the overload resolutions
515 static std::unordered_map<TClass *, Func> dispatchMap;
516
517 auto found = dispatchMap.find(tclass);
518
519 if (found != dispatchMap.end()) {
520 func = found->second;
521 } else {
522 // Can probably done with CppInterop in the future to avoid string manipulation.
523 std::stringstream cmd;
524 cmd << "&RooFit::Experimental::CodegenIntegralImplCaller<" << tclass->GetName() << ">::call;";
525 func = reinterpret_cast<Func>(gInterpreter->ProcessLine(cmd.str().c_str()));
526 dispatchMap[tclass] = func;
527 }
528
529 return func(arg, code, rangeName, ctx);
530}
531
532} // namespace
533
535{
536 if (arg.numIntCatVars().empty() && arg.numIntRealVars().empty()) {
537 ctx.addResult(&arg, codegenIntegral(const_cast<RooAbsReal &>(arg.integrand()), arg.mode(), arg.intRange(), ctx));
538 return;
539 }
540
541 if (arg.intVars().size() != 1 || arg.numIntRealVars().size() != 1) {
542 std::stringstream errorMsg;
543 errorMsg << "Only analytical integrals and 1D numeric integrals are supported for AD for class"
544 << arg.integrand().GetName();
545 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
546 throw std::runtime_error(errorMsg.str().c_str());
547 }
548
549 auto &intVar = static_cast<RooAbsRealLValue &>(*arg.numIntRealVars()[0]);
550
551 std::string obsName = ctx.getTmpVarName();
552 std::string oldIntVarResult = ctx.getResult(intVar);
553 ctx.addResult(&intVar, "obs[0]");
554
555 std::string funcName = ctx.buildFunction(arg.integrand(), {});
556
557 std::stringstream ss;
558
559 ss << "double " << obsName << "[1];\n";
560
561 std::string resName = RooFit::Detail::makeValidVarName(arg.GetName()) + "Result";
562 ctx.addResult(&arg, resName);
563 ctx.addToGlobalScope("double " + resName + " = 0.0;\n");
564
565 // TODO: once Clad has support for higher-order functions (follow also the
566 // Clad issue #637), we could refactor this code into an actual function
567 // instead of hardcoding it here as a string.
568 ss << "{\n"
569 << " const int n = 1000; // number of sampling points\n"
570 << " double d = " << intVar.getMax(arg.intRange()) << " - " << intVar.getMin(arg.intRange()) << ";\n"
571 << " double eps = d / n;\n"
572 << " for (int i = 0; i < n; ++i) {\n"
573 << " " << obsName << "[0] = " << intVar.getMin(arg.intRange()) << " + eps * i;\n"
574 << " double tmpA = " << funcName << "(params, " << obsName << ", xlArr);\n"
575 << " " << obsName << "[0] = " << intVar.getMin(arg.intRange()) << " + eps * (i + 1);\n"
576 << " double tmpB = " << funcName << "(params, " << obsName << ", xlArr);\n"
577 << " " << resName << " += (tmpA + tmpB) * 0.5 * eps;\n"
578 << " }\n"
579 << "}\n";
580
581 ctx.addToGlobalScope(ss.str());
582
583 ctx.addResult(&intVar, oldIntVarResult);
584}
585
587{
588 ctx.addResult(&arg, realSumPdfTranslateImpl(ctx, arg, arg.funcList(), arg.coefList(), false));
589}
590
592{
593 ctx.addResult(&arg, realSumPdfTranslateImpl(ctx, arg, arg.funcList(), arg.coefList(), false));
594}
595
597{
598 if (!arg.isConstant()) {
599 ctx.addResult(&arg, arg.GetName());
600 }
601 // Just return a formatted version of the const value.
602 // Formats to the maximum precision.
603 constexpr auto max_precision{std::numeric_limits<double>::digits10 + 1};
604 std::stringstream ss;
605 ss.precision(max_precision);
606 // Just use toString to make sure we do not output 'inf'.
607 // This is really ugly for large numbers...
608 ss << std::fixed << RooNumber::toString(arg.getVal());
609 ctx.addResult(&arg, ss.str());
610}
611
613{
614 ctx.addResult(&arg, ctx.buildCall(mathFunc("recursiveFraction"), arg.variables(), arg.variables().size()));
615}
616
618{
619 auto const &interpCodes = arg.interpolationCodes();
620
621 unsigned int n = interpCodes.size();
622
623 int interpCode = interpCodes[0];
624 // To get consistent codes with the PiecewiseInterpolation
625 if (interpCode == 4) {
626 interpCode = 5;
627 }
628
629 for (unsigned int i = 1; i < n; i++) {
630 if (interpCodes[i] != interpCodes[0]) {
632 << "FlexibleInterpVar::evaluate ERROR: Code Squashing AD does not yet support having "
633 "different interpolation codes for the same class object "
634 << std::endl;
635 }
636 }
637
638 std::string const &resName = ctx.buildCall(mathFunc("flexibleInterp"), interpCode, arg.variables(), n, arg.low(),
639 arg.high(), arg.globalBoundary(), arg.nominal(), 1.0);
640 ctx.addResult(&arg, resName);
641}
642
644{
645 ctx.addResult(&arg, "1.0");
646}
647
648////////////////////////////////////////////////////////////////////////////////
649/// This function defines the analytical integral translation for the class.
650///
651/// \param[in] code The code that decides the integrands.
652/// \param[in] rangeName Name of the normalization range.
653/// \param[in] ctx An object to manage auxiliary information for code-squashing.
654///
655/// \returns The representative code string of the integral for the given object.
656std::string codegenIntegralImpl(RooAbsReal &arg, int, const char *, CodegenContext &)
657{
658 std::stringstream errorMsg;
659 errorMsg << "An analytical integral function for class \"" << arg.ClassName() << "\" has not yet been implemented.";
660 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
661 throw std::runtime_error(errorMsg.str().c_str());
662}
663
664std::string codegenIntegralImpl(RooBernstein &arg, int, const char *rangeName, CodegenContext &ctx)
665{
666 arg.fillBuffer(); // to get the right xmin() and xmax()
667 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
668 return ctx.buildCall(mathFunc("bernsteinIntegral"), x.getMin(rangeName), x.getMax(rangeName), arg.xmin(), arg.xmax(),
669 arg.coefList(), arg.coefList().size());
670}
671
672std::string codegenIntegralImpl(RooBifurGauss &arg, int code, const char *rangeName, CodegenContext &ctx)
673{
674 auto &constant = code == 1 ? arg.getMean() : arg.getX();
675 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(code == 1 ? arg.getX() : arg.getMean());
676
677 return ctx.buildCall(mathFunc("bifurGaussIntegral"), integrand.getMin(rangeName), integrand.getMax(rangeName),
678 constant, arg.getSigmaL(), arg.getSigmaR());
679}
680
681std::string codegenIntegralImpl(RooCBShape &arg, int /*code*/, const char *rangeName, CodegenContext &ctx)
682{
683 auto &m = dynamic_cast<RooAbsRealLValue const &>(arg.getM());
684 return ctx.buildCall(mathFunc("cbShapeIntegral"), m.getMin(rangeName), m.getMax(rangeName), arg.getM0(),
685 arg.getSigma(), arg.getAlpha(), arg.getN());
686}
687
688std::string codegenIntegralImpl(RooChebychev &arg, int, const char *rangeName, CodegenContext &ctx)
689{
690 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
691 double xmax = x.getMax(arg.refRangeName());
692 double xmin = x.getMin(arg.refRangeName());
693 unsigned int sz = arg.coefList().size();
694
695 return ctx.buildCall(mathFunc("chebychevIntegral"), arg.coefList(), sz, xmin, xmax, x.getMin(rangeName),
696 x.getMax(rangeName));
697}
698
699std::string codegenIntegralImpl(RooEfficiency &, int, const char *, CodegenContext &)
700{
701 return "1.0";
702}
703
704std::string codegenIntegralImpl(RooExponential &arg, int code, const char *rangeName, CodegenContext &ctx)
705{
706 bool isOverX = code == 1;
707
708 std::string constant;
709 if (arg.negateCoefficient() && isOverX) {
710 constant += "-";
711 }
712 constant += ctx.getResult(isOverX ? arg.coefficient() : arg.variable());
713
714 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(isOverX ? arg.variable() : arg.coefficient());
715
716 double min = integrand.getMin(rangeName);
717 double max = integrand.getMax(rangeName);
718
719 if (!isOverX && arg.negateCoefficient()) {
720 std::swap(min, max);
721 min = -min;
722 max = -max;
723 }
724
725 return ctx.buildCall(mathFunc("exponentialIntegral"), min, max, constant);
726}
727
728std::string codegenIntegralImpl(RooGamma &arg, int, const char *rangeName, CodegenContext &ctx)
729{
730 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.getX());
731 const std::string a =
732 ctx.buildCall("ROOT::Math::gamma_cdf", x.getMax(rangeName), arg.getGamma(), arg.getBeta(), arg.getMu());
733 const std::string b =
734 ctx.buildCall("ROOT::Math::gamma_cdf", x.getMin(rangeName), arg.getGamma(), arg.getBeta(), arg.getMu());
735 return a + " - " + b;
736}
737
738std::string codegenIntegralImpl(RooGaussian &arg, int code, const char *rangeName, CodegenContext &ctx)
739{
740 auto &constant = code == 1 ? arg.getMean() : arg.getX();
741 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(code == 1 ? arg.getX() : arg.getMean());
742
743 return ctx.buildCall(mathFunc("gaussianIntegral"), integrand.getMin(rangeName), integrand.getMax(rangeName),
744 constant, arg.getSigma());
745}
746
747namespace {
748
749std::string rooHistIntegralTranslateImpl(int code, RooAbsArg const &arg, RooDataHist const &dataHist,
750 const RooArgSet &obs, bool histFuncMode)
751{
752 if (((2 << obs.size()) - 1) != code) {
753 oocoutE(&arg, InputArguments) << "RooHistPdf::integral(" << arg.GetName()
754 << ") ERROR: AD currently only supports integrating over all histogram observables."
755 << std::endl;
756 return "";
757 }
758 return std::to_string(dataHist.sum(histFuncMode));
759}
760
761} // namespace
762
763std::string codegenIntegralImpl(RooHistFunc &arg, int code, const char *, CodegenContext &)
764{
765 return rooHistIntegralTranslateImpl(code, arg, arg.dataHist(), arg.variables(), true);
766}
767
768std::string codegenIntegralImpl(RooHistPdf &arg, int code, const char *, CodegenContext &)
769{
770 return rooHistIntegralTranslateImpl(code, arg, arg.dataHist(), arg.variables(), false);
771}
772
773std::string codegenIntegralImpl(RooLandau &arg, int, const char *rangeName, CodegenContext &ctx)
774{
775 // Don't do anything with "code". It can only be "1" anyway (see
776 // implementation of getAnalyticalIntegral).
777 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.getX());
778 const std::string a = ctx.buildCall("ROOT::Math::landau_cdf", x.getMax(rangeName), arg.getSigma(), arg.getMean());
779 const std::string b = ctx.buildCall("ROOT::Math::landau_cdf", x.getMin(rangeName), arg.getSigma(), arg.getMean());
780 return ctx.getResult(arg.getSigma()) + " * " + "(" + a + " - " + b + ")";
781}
782
783std::string codegenIntegralImpl(RooLognormal &arg, int, const char *rangeName, CodegenContext &ctx)
784{
785 std::string funcName = arg.useStandardParametrization() ? "logNormalIntegralStandard" : "logNormalIntegral";
786 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.getX());
787 return ctx.buildCall(mathFunc(funcName), x.getMin(rangeName), x.getMax(rangeName), arg.getMedian(), arg.getShapeK());
788}
789
790std::string codegenIntegralImpl(RooMultiVarGaussian &arg, int code, const char *rangeName, CodegenContext &)
791{
792 if (code != -1) {
793 std::stringstream errorMsg;
794 errorMsg << "Partial integrals over RooMultiVarGaussian are not supported.";
795 oocoutE(&arg, Minimization) << errorMsg.str() << std::endl;
796 throw std::runtime_error(errorMsg.str().c_str());
797 }
798
799 return std::to_string(arg.analyticalIntegral(code, rangeName));
800}
801
802std::string codegenIntegralImpl(RooPoisson &arg, int code, const char *rangeName, CodegenContext &ctx)
803{
804 assert(code == 1 || code == 2);
805 std::string xName = ctx.getResult(arg.getX());
806 if (!arg.getNoRounding())
807 xName = "std::floor(" + xName + ")";
808
809 auto &integrand = dynamic_cast<RooAbsRealLValue const &>(code == 1 ? arg.getX() : arg.getMean());
810 // Since the integral function is the same for both codes, we need to make sure the indexed observables do not appear
811 // in the function if they are not required.
812 xName = code == 1 ? "0" : xName;
813 return ctx.buildCall(mathFunc("poissonIntegral"), code, arg.getMean(), xName, integrand.getMin(rangeName),
814 integrand.getMax(rangeName), arg.getProtectNegativeMean());
815}
816
817std::string codegenIntegralImpl(RooPolyVar &arg, int, const char *rangeName, CodegenContext &ctx)
818{
819 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
820 const double xmin = x.getMin(rangeName);
821 const double xmax = x.getMax(rangeName);
822 const unsigned sz = arg.coefList().size();
823 if (!sz)
824 return std::to_string(arg.lowestOrder() ? xmax - xmin : 0.0);
825
826 return ctx.buildCall(mathFunc("polynomialIntegral"), arg.coefList(), sz, arg.lowestOrder(), xmin, xmax);
827}
828
829std::string codegenIntegralImpl(RooPolynomial &arg, int, const char *rangeName, CodegenContext &ctx)
830{
831 auto &x = dynamic_cast<RooAbsRealLValue const &>(arg.x());
832 const double xmin = x.getMin(rangeName);
833 const double xmax = x.getMax(rangeName);
834 const unsigned sz = arg.coefList().size();
835 if (!sz)
836 return std::to_string(arg.lowestOrder() ? xmax - xmin : 0.0);
837
838 return ctx.buildCall(mathFunc("polynomialIntegral<true>"), arg.coefList(), sz, arg.lowestOrder(), xmin, xmax);
839}
840
841std::string codegenIntegralImpl(RooRealSumPdf &arg, int code, const char *rangeName, CodegenContext &ctx)
842{
843 // Re-use translate, since integration is linear.
844 return realSumPdfTranslateImpl(ctx, arg, arg.funcIntListFromCache(code, rangeName), arg.coefList(), false);
845}
846
847std::string codegenIntegralImpl(RooUniform &arg, int code, const char *rangeName, CodegenContext &)
848{
849 // The integral of a uniform distribution is static, so we can just hardcode
850 // the result in a string.
851 return std::to_string(arg.analyticalIntegral(code, rangeName));
852}
853
854} // namespace Experimental
855} // namespace RooFit
#define b(i)
Definition RSha256.hxx:100
#define a(i)
Definition RSha256.hxx:99
#define oocoutE(o, a)
#define ooccoutE(o, a)
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
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
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
bool isConstant() const
Check if the "Constant" attribute is set.
Definition RooAbsArg.h:304
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
virtual Int_t numEntries() const
Return number of entries in dataset, i.e., count unweighted entries.
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:59
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:103
TClass * IsA() const override
Definition RooAbsReal.h:547
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
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(RooAbsArg const *klass, RooFit::Experimental::CodegenContext &ctx, const RooAbsCollection &coords, bool reverse=false) const
double binVolume(std::size_t i) const
Return bin volume of i-th bin.
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:206
RooProdPdf::CacheElem const & cache() const
Definition RooProdPdf.h:250
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.
std::string buildFunction(RooAbsArg const &arg, std::map< RooFit::Detail::DataKey, std::size_t > const &outputSizes={})
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.
std::string buildArg(RooAbsCollection const &x)
Function to save a RooListProxy as an array in the squashed code.
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
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 propability 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.
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
static std::string toString(double x)
Returns an std::to_string compatible number (i.e.
Definition RooNumber.cxx:31
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
RooRealProxy const & x() const
Definition RooPolyVar.h:37
RooArgList const & coefList() const
Definition RooPolyVar.h:38
int lowestOrder() const
Definition RooPolyVar.h:39
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.
std::unique_ptr< RooAbsReal > _rearrangedNum
Definition RooProdPdf.h:118
std::unique_ptr< RooAbsReal > _rearrangedDen
Definition RooProdPdf.h:119
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.
TClass instances represent classes, structs and namespaces in the ROOT type system.
Definition TClass.h:81
const Element * GetMatrixArray() const override
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
virtual const char * ClassName() const
Returns name of class to which the object belongs.
Definition TObject.cxx:225
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)
std::string codegenIntegralImpl(RooAbsReal &arg, int code, const char *rangeName, CodegenContext &ctx)
This function defines the analytical integral translation for the class.
The namespace RooFit contains mostly switches that change the behaviour of functions of PDFs (or othe...
Definition CodegenImpl.h:64
@ InputArguments
TMarker m
Definition textangle.C:8
static uint64_t sum(uint64_t i)
Definition Factory.cxx:2345