Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooLegacyExpPoly.cxx
Go to the documentation of this file.
1/*
2 * Project: RooFit
3 *
4 * Copyright (c) 2022, CERN
5 *
6 * Redistribution and use in source and binary forms,
7 * with or without modification, are permitted according to the terms
8 * listed in LICENSE (http://roofit.sourceforge.net/license.txt)
9 */
10
11/** \class RooLegacyExpPoly
12 \ingroup Roofit
13
14RooLegacyExpPoly implements a polynomial PDF of the form \f[ f(x) =
15\mathcal{N} \cdot \exp( \sum_{i} a_{i} * x^{i} ) \f] \f$ \mathcal{N}
16\f$ is a normalisation constant that is automatically calculated when
17the function is used in computations.
18
19The sum can be truncated at the low end. See the main constructor
20RooLegacyExpPoly::RooLegacyExpPoly(const char*, const char*, RooAbsReal&, const RooArgList&, int)
21
22\image html RooLegacyExpPoly.png
23
24**/
25
26#include <RooLegacyExpPoly.h>
27
28#include <RooAbsReal.h>
29#include <RooArgList.h>
30#include <RooMath.h>
31#include <RooMsgService.h>
32#include <RooRealVar.h>
33#include "RooBatchCompute.h"
34
35#include <TMath.h>
36#include <TError.h>
37
38#include <array>
39#include <cassert>
40#include <cmath>
41#include <complex>
42#include <sstream>
43
45
46////////////////////////////////////////////////////////////////////////////////
47/// Create a polynomial in the variable `x`.
48/// \param[in] name Name of the PDF
49/// \param[in] title Title for plotting the PDF
50/// \param[in] x The variable of the polynomial
51/// \param[in] coefList The coefficients \f$ a_i \f$
52/// \param[in] lowestOrder [optional] Truncate the sum such that it skips the lower orders:
53/// \f[
54/// 1. + \sum_{i=0}^{\mathrm{coefList.size()}} a_{i} * x^{(i + \mathrm{lowestOrder})}
55/// \f]
56///
57/// This means that
58/// \code{.cpp}
59/// RooLegacyExpPoly pol("pol", "pol", x, RooArgList(a, b), lowestOrder = 2)
60/// \endcode
61/// computes
62/// \f[
63/// \mathrm{pol}(x) = 1 * x^0 + (0 * x^{\ldots}) + a * x^2 + b * x^3.
64/// \f]
65
66RooLegacyExpPoly::RooLegacyExpPoly(const char *name, const char *title, RooAbsReal &x, const RooArgList &coefList, int lowestOrder)
67 : RooAbsPdf(name, title),
68 _x("x", "Dependent", this, x),
69 _coefList("coefList", "List of coefficients", this),
70 _lowestOrder(lowestOrder)
71{
72 // Check lowest order
73 if (_lowestOrder < 0) {
74 coutE(InputArguments) << "RooLegacyExpPoly::ctor(" << GetName()
75 << ") WARNING: lowestOrder must be >=0, setting value to 0" << std::endl;
76 _lowestOrder = 0;
77 }
78
80}
81
82////////////////////////////////////////////////////////////////////////////////
83/// Copy constructor
84
86 : RooAbsPdf(other, name),
87 _x("x", this, other._x),
88 _coefList("coefList", this, other._coefList),
89 _lowestOrder(other._lowestOrder)
90{
91}
92
93////////////////////////////////////////////////////////////////////////////////
94
96{
97 // Calculate and return value of polynomial
98
99 const unsigned sz = _coefList.size();
100 const int lowestOrder = _lowestOrder;
101 if (!sz)
102 return lowestOrder ? 1. : 0.;
103 std::vector<double> coefs;
104 coefs.reserve(sz);
105
106 const RooArgSet *nset = _coefList.nset();
107 for (auto coef : _coefList) {
108 coefs.push_back(static_cast<RooAbsReal *>(coef)->getVal(nset));
109 };
110 const double x = _x;
111 double xpow = std::pow(x, lowestOrder);
112 double retval = 0;
113 for (size_t i = 0; i < sz; ++i) {
114 retval += coefs[i] * xpow;
115 xpow *= x;
116 }
117
118 if (std::numeric_limits<double>::max_exponent < retval) {
119 coutE(InputArguments) << "RooLegacyExpPoly::evaluateLog(" << GetName() << ") ERROR: exponent at " << x
120 << " larger than allowed maximum, result will be infinite! " << retval << " > "
121 << std::numeric_limits<double>::max_exponent << " in " << this->getFormulaExpression(true)
122 << std::endl;
123 }
124 return retval;
125}
126
127////////////////////////////////////////////////////////////////////////////////
128
129/// Compute multiple values of ExpPoly distribution.
131{
132 std::vector<std::span<const double>> vars;
133 vars.reserve(_coefList.size() + 1);
134 vars.push_back(ctx.at(_x));
135
136 std::vector<double> coefVals;
137 for (RooAbsArg *coef : _coefList) {
138 vars.push_back(ctx.at(coef));
139 }
140
141 std::array<double, 2> args{static_cast<double>(_lowestOrder), static_cast<double>(_coefList.size())};
142
143 RooBatchCompute::compute(ctx.config(this), RooBatchCompute::ExpPoly, ctx.output(), vars, args);
144}
145
146////////////////////////////////////////////////////////////////////////////////
147
149{
150 // Adjust the limits of all the coefficients to reflect the numeric boundaries
151
152 const unsigned sz = _coefList.size();
153 double max = std::numeric_limits<double>::max_exponent / sz;
154 const int lowestOrder = _lowestOrder;
155 std::vector<double> coefs;
156 coefs.reserve(sz);
157
158 RooRealVar *x = dynamic_cast<RooRealVar *>(&(*_x));
159 if (x) {
160 const double xmax = x->getMax();
161 double xmaxpow = std::pow(xmax, lowestOrder);
162 for (size_t i = 0; i < sz; ++i) {
163 double thismax = max / xmaxpow;
164 RooRealVar *coef = dynamic_cast<RooRealVar *>(this->_coefList.at(i));
165 if (coef) {
166 coef->setVal(thismax);
167 coef->setMax(thismax);
168 }
169 xmaxpow *= xmax;
170 }
171 }
172}
173
174////////////////////////////////////////////////////////////////////////////////
175
177{
178 // Calculate and return value of function
179
180 const double logval = this->evaluateLog();
181 const double val = std::exp(logval);
182 if (std::isinf(val)) {
183 coutE(InputArguments) << "RooLegacyExpPoly::evaluate(" << GetName()
184 << ") ERROR: result of exponentiation is infinite! exponent was " << logval << std::endl;
185 }
186 return val;
187}
188
189////////////////////////////////////////////////////////////////////////////////
190
191double RooLegacyExpPoly::getLogVal(const RooArgSet *nset) const
192{
193 return RooAbsPdf::getLogVal(nset);
194 // return this->evaluateLog();
195}
196
197////////////////////////////////////////////////////////////////////////////////
198
199int RooLegacyExpPoly::getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char * /*rangeName*/) const
200{
201
202 if ((_coefList.size() + _lowestOrder < 4) &&
203 ((_coefList.size() + _lowestOrder < 3) ||
204 (static_cast<RooAbsReal *>(_coefList.at(2 - _lowestOrder))->getVal() <= 0)) &&
205 matchArgs(allVars, analVars, _x)) {
206 return 0;
207 }
208 return 0;
209}
210
211////////////////////////////////////////////////////////////////////////////////
212
213#define PI TMath::Pi()
214
215namespace {
216double deltaerf(long double x1, long double x2)
217{
218 // several things happening here
219 // 1. use "erf(x) = -erf(-x)" to evaluate the function only ever for the positive side (higher precision)
220 // 2. use "erf(x) = 1.-erfc(x)" and, instead of "erf(x1) - erf(x2)", do "(1.-erfc(x1)) - (1.-erfc(x2)) = erfc(x2) -
221 // erfc(x1)"
222 double y2 = x1 > 0 ? erfc(x1) : -erfc(-x1);
223 double y1 = x2 > 0 ? erfc(x2) : -erfc(-x2);
224 if (y1 == y2) {
225 std::cout << "WARNING in calculation of analytical integral limited by numerical precision" << std::endl;
226 std::cout << "x: " << x1 << " , " << x2 << std::endl;
227 std::cout << "y: " << y1 << " , " << y2 << std::endl;
228 }
229 return y1 - y2;
230}
231
232double deltaerfi(double x1, double x2)
233{
234 std::complex<double> u1 = {x1, 0.};
235 std::complex<double> u2 = {x2, 0.};
236
237 std::complex<double> y2 = x2 > 0 ? RooMath::faddeeva(u2) : RooMath::faddeeva(-u2);
238 std::complex<double> y1 = x1 > 0 ? RooMath::faddeeva(u1) : RooMath::faddeeva(-u1);
239 if (y1 == y2) {
240 std::cout << "WARNING in calculation of analytical integral limited by numerical precision" << std::endl;
241 std::cout << "x: " << x1 << " , " << x2 << std::endl;
242 std::cout << "y: " << y1 << " , " << y2 << std::endl;
243 }
244 return y1.imag() - y2.imag();
245}
246} // namespace
247
248double RooLegacyExpPoly::analyticalIntegral(int /*code*/, const char *rangeName) const
249{
250 const double xmin = _x.min(rangeName);
251 const double xmax = _x.max(rangeName);
252 const unsigned sz = _coefList.size();
253 if (!sz)
254 return xmax - xmin;
255
256 std::vector<double> coefs;
257 coefs.reserve(sz);
258 const RooArgSet *nset = _coefList.nset();
259 for (auto c : _coefList) {
260 coefs.push_back(static_cast<RooAbsReal *>(c)->getVal(nset));
261 }
262
263 switch (_coefList.size() + _lowestOrder) {
264 case 1: return xmax - xmin;
265 case 2: {
266 const double a = coefs[1 - _lowestOrder];
267 if (a != 0) {
268 return 1. / a * (exp(a * xmax) - exp(a * xmin)) * (_lowestOrder == 0 ? exp(coefs[0]) : 1);
269 } else {
270 return xmax - xmin;
271 }
272 }
273 case 3: {
274 const double a = coefs[2 - _lowestOrder];
275 const double b = _lowestOrder == 2 ? 0. : coefs[1 - _lowestOrder];
276 const double c = _lowestOrder == 0 ? coefs[0] : 0.;
277 const double absa = std::abs(a);
278 const double sqrta = std::sqrt(absa);
279 if (a < 0) {
280 double d = ::deltaerf((-b + 2 * absa * xmax) / (2 * sqrta), (-b + 2 * absa * xmin) / (2 * sqrta));
281 double retval = exp(b * b / (4 * absa) + c) * std::sqrt(PI) * d / (2 * sqrta);
282 return retval;
283 } else if (a > 0) {
284 double d = ::deltaerfi((b + 2 * absa * xmax) / (2 * sqrta), (b + 2 * absa * xmin) / (2 * sqrta));
285 double retval = exp(-b * b / (4 * absa) + c) * std::sqrt(PI) * d / (2 * sqrta);
286 return retval;
287 } else if (b != 0) {
288 return 1. / b * (std::exp(b * xmax) - exp(b * xmin)) * exp(c);
289 } else {
290 return xmax - xmin;
291 }
292 }
293 }
294 return 0.;
295}
296
297////////////////////////////////////////////////////////////////////////////////
298
299std::string RooLegacyExpPoly::getFormulaExpression(bool expand) const
300{
301 std::stringstream ss;
302 ss << "exp(";
303 int order = _lowestOrder;
304 for (auto coef : _coefList) {
305 if (order != _lowestOrder)
306 ss << "+";
307 if (expand) {
308 ss << (static_cast<RooAbsReal *>(coef))->getVal();
309 } else {
310 ss << coef->GetName();
311 }
312 ss << "*pow(" << _x.GetName() << "," << order << ")";
313 ++order;
314 }
315 ss << ")";
316 return ss.str();
317}
#define d(i)
Definition RSha256.hxx:102
#define b(i)
Definition RSha256.hxx:100
#define c(i)
Definition RSha256.hxx:101
#define a(i)
Definition RSha256.hxx:99
#define coutE(a)
#define ClassImp(name)
Definition Rtypes.h:382
#define PI
Option_t Option_t TPoint TPoint const char x2
Option_t Option_t TPoint TPoint const char x1
Option_t Option_t TPoint TPoint const char y2
Option_t Option_t TPoint TPoint const char y1
char name[80]
Definition TGX11.cxx:110
float xmin
float xmax
Common abstract base class for objects that represent a value and a "shape" in RooFit.
Definition RooAbsArg.h:77
Storage_t::size_type size() const
bool addTyped(const RooAbsCollection &list, bool silent=false)
Adds elements of a given RooAbsCollection to the container if they match the specified type.
Abstract interface for all probability density functions.
Definition RooAbsPdf.h:40
virtual double getLogVal(const RooArgSet *set=nullptr) const
Return the log of the current value with given normalization An error message is printed if the argum...
const RooArgSet * nset() const
Definition RooAbsProxy.h:52
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
bool matchArgs(const RooArgSet &allDeps, RooArgSet &numDeps, const RooArgProxy &a) const
Utility function for use in getAnalyticalIntegral().
RooArgList is a container object that can hold multiple RooAbsArg objects.
Definition RooArgList.h:22
RooAbsArg * at(Int_t idx) const
Return object at given index, or nullptr if index is out of range.
Definition RooArgList.h:110
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:24
std::span< const double > at(RooAbsArg const *arg, RooAbsArg const *caller=nullptr)
std::span< double > output()
RooBatchCompute::Config config(RooAbsArg const *arg) const
RooLegacyExpPoly implements a polynomial PDF of the form.
int getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
double evaluateLog() const
RooArgList const & coefList() const
Get the coefficient list.
RooListProxy _coefList
RooAbsReal const & x() const
Get the x variable.
double getLogVal(const RooArgSet *nset) const override
Return the log of the current value with given normalization An error message is printed if the argum...
std::string getFormulaExpression(bool expand) const
double evaluate() const override
Evaluation.
void doEval(RooFit::EvalContext &) const override
Compute multiple values of ExpPoly distribution.
double analyticalIntegral(int code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
int lowestOrder() const
Return the order for the first coefficient in the list.
static std::complex< double > faddeeva(std::complex< double > z)
evaluate Faddeeva function for complex argument
Definition RooMath.cxx:30
Variable that can be changed from the outside.
Definition RooRealVar.h:37
void setVal(double value) override
Set value of variable to 'value'.
void setMax(const char *name, double value)
Set maximum of name range to given value.
double max(const char *rname=nullptr) const
Query upper limit of range. This requires the payload to be RooAbsRealLValue or derived.
double min(const char *rname=nullptr) const
Query lower limit of range. This requires the payload to be RooAbsRealLValue or derived.
const char * GetName() const override
Returns name of object.
Definition TNamed.h:47
Double_t x[n]
Definition legend1.C:17
void compute(Config cfg, Computer comp, std::span< double > output, VarSpan vars, ArgSpan extraArgs={})