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