Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooExpPoly.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 RooExpPoly
12 \ingroup Roofit
13
14RooExpPoly 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
20RooExpPoly::RooExpPoly(const char*, const char*, RooAbsReal&, const RooArgList&, int)
21
22\image html RooExpPoly.png
23
24**/
25
26#include <RooExpPoly.h>
27
28#include <RooAbsReal.h>
29#include <RooArgList.h>
30#include <RooMath.h>
31#include <RooMsgService.h>
32#include <RooRealVar.h>
33
34#include <TMath.h>
35#include <TError.h>
36
37#include <cmath>
38#include <sstream>
39#include <cassert>
40#include <complex>
41
43
44////////////////////////////////////////////////////////////////////////////////
45/// Create a polynomial in the variable `x`.
46/// \param[in] name Name of the PDF
47/// \param[in] title Title for plotting the PDF
48/// \param[in] x The variable of the polynomial
49/// \param[in] coefList The coefficients \f$ a_i \f$
50/// \param[in] lowestOrder [optional] Truncate the sum such that it skips the lower orders:
51/// \f[
52/// 1. + \sum_{i=0}^{\mathrm{coefList.size()}} a_{i} * x^{(i + \mathrm{lowestOrder})}
53/// \f]
54///
55/// This means that
56/// \code{.cpp}
57/// RooExpPoly pol("pol", "pol", x, RooArgList(a, b), lowestOrder = 2)
58/// \endcode
59/// computes
60/// \f[
61/// \mathrm{pol}(x) = 1 * x^0 + (0 * x^{\ldots}) + a * x^2 + b * x^3.
62/// \f]
63
64RooExpPoly::RooExpPoly(const char *name, const char *title, RooAbsReal &x, const RooArgList &coefList, int lowestOrder)
65 : RooAbsPdf(name, title),
66 _x("x", "Dependent", this, x),
67 _coefList("coefList", "List of coefficients", this),
68 _lowestOrder(lowestOrder)
69{
70 // Check lowest order
71 if (_lowestOrder < 0) {
72 coutE(InputArguments) << "RooExpPoly::ctor(" << GetName()
73 << ") WARNING: lowestOrder must be >=0, setting value to 0" << std::endl;
74 _lowestOrder = 0;
75 }
76
77 for (auto coef : coefList) {
78 if (!dynamic_cast<RooAbsReal *>(coef)) {
79 coutE(InputArguments) << "RooExpPoly::ctor(" << GetName() << ") ERROR: coefficient " << coef->GetName()
80 << " is not of type RooAbsReal" << std::endl;
81 R__ASSERT(0);
82 }
83 _coefList.add(*coef);
84 }
85}
86
87////////////////////////////////////////////////////////////////////////////////
88/// Copy constructor
89
90RooExpPoly::RooExpPoly(const RooExpPoly &other, const char *name)
91 : RooAbsPdf(other, name),
92 _x("x", this, other._x),
93 _coefList("coefList", this, other._coefList),
94 _lowestOrder(other._lowestOrder)
95{
96}
97
98////////////////////////////////////////////////////////////////////////////////
99
101{
102 // Calculate and return value of polynomial
103
104 const unsigned sz = _coefList.size();
105 const int lowestOrder = _lowestOrder;
106 if (!sz)
107 return lowestOrder ? 1. : 0.;
108 std::vector<double> coefs;
109 coefs.reserve(sz);
110
111 const RooArgSet *nset = _coefList.nset();
112 for (auto coef : _coefList) {
113 coefs.push_back(static_cast<RooAbsReal *>(coef)->getVal(nset));
114 };
115 const double x = _x;
116 double xpow = std::pow(x, lowestOrder);
117 double retval = 0;
118 for (size_t i = 0; i < sz; ++i) {
119 retval += coefs[i] * xpow;
120 xpow *= x;
121 }
122
123 if (std::numeric_limits<double>::max_exponent < retval) {
124 coutE(InputArguments) << "RooExpPoly::evaluateLog(" << GetName() << ") ERROR: exponent at " << x
125 << " larger than allowed maximum, result will be infinite! " << retval << " > "
126 << std::numeric_limits<double>::max_exponent << " in " << this->getFormulaExpression(true)
127 << std::endl;
128 }
129 return retval;
130}
131
132////////////////////////////////////////////////////////////////////////////////
133
135{
136 // Adjust the limits of all the coefficients to reflect the numeric boundaries
137
138 const unsigned sz = _coefList.size();
139 double max = std::numeric_limits<double>::max_exponent / sz;
140 const int lowestOrder = _lowestOrder;
141 std::vector<double> coefs;
142 coefs.reserve(sz);
143
144 RooRealVar *x = dynamic_cast<RooRealVar *>(&(*_x));
145 if (x) {
146 const double xmax = x->getMax();
147 double xmaxpow = std::pow(xmax, lowestOrder);
148 for (size_t i = 0; i < sz; ++i) {
149 double thismax = max / xmaxpow;
150 RooRealVar *coef = dynamic_cast<RooRealVar *>(this->_coefList.at(i));
151 if (coef) {
152 coef->setVal(thismax);
153 coef->setMax(thismax);
154 }
155 xmaxpow *= xmax;
156 }
157 }
158}
159
160////////////////////////////////////////////////////////////////////////////////
161
163{
164 // Calculate and return value of function
165
166 const double logval = this->evaluateLog();
167 const double val = exp(logval);
168 if (std::isinf(val)) {
169 coutE(InputArguments) << "RooExpPoly::evaluate(" << GetName()
170 << ") ERROR: result of exponentiation is infinite! exponent was " << logval << std::endl;
171 }
172 return val;
173}
174
175////////////////////////////////////////////////////////////////////////////////
176
177double RooExpPoly::getLogVal(const RooArgSet *nset) const
178{
179 return RooAbsPdf::getLogVal(nset);
180 // return this->evaluateLog();
181}
182
183////////////////////////////////////////////////////////////////////////////////
184
185int RooExpPoly::getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char * /*rangeName*/) const
186{
187
188 if ((_coefList.size() + _lowestOrder < 4) &&
189 ((_coefList.size() + _lowestOrder < 3) ||
190 (static_cast<RooAbsReal *>(_coefList.at(2 - _lowestOrder))->getVal() <= 0)) &&
191 matchArgs(allVars, analVars, _x)) {
192 return 0;
193 }
194 return 0;
195}
196
197////////////////////////////////////////////////////////////////////////////////
198
199#define PI TMath::Pi()
200
201namespace {
202double deltaerf(long double x1, long double x2)
203{
204 // several things happening here
205 // 1. use "erf(x) = -erf(-x)" to evaluate the function only ever for the positive side (higher precision)
206 // 2. use "erf(x) = 1.-erfc(x)" and, instead of "erf(x1) - erf(x2)", do "(1.-erfc(x1)) - (1.-erfc(x2)) = erfc(x2) -
207 // erfc(x1)"
208 double y2 = x1 > 0 ? erfc(x1) : -erfc(-x1);
209 double y1 = x2 > 0 ? erfc(x2) : -erfc(-x2);
210 if (y1 == y2) {
211 std::cout << "WARNING in calculation of analytical integral limited by numerical precision" << std::endl;
212 std::cout << "x: " << x1 << " , " << x2 << std::endl;
213 std::cout << "y: " << y1 << " , " << y2 << std::endl;
214 }
215 return y1 - y2;
216}
217
218double deltaerfi(double x1, double x2)
219{
220 std::complex<double> u1 = {x1, 0.};
221 std::complex<double> u2 = {x2, 0.};
222
223 std::complex<double> y2 = x2 > 0 ? RooMath::faddeeva(u2) : RooMath::faddeeva(-u2);
224 std::complex<double> y1 = x1 > 0 ? RooMath::faddeeva(u1) : RooMath::faddeeva(-u1);
225 if (y1 == y2) {
226 std::cout << "WARNING in calculation of analytical integral limited by numerical precision" << std::endl;
227 std::cout << "x: " << x1 << " , " << x2 << std::endl;
228 std::cout << "y: " << y1 << " , " << y2 << std::endl;
229 }
230 return y1.imag() - y2.imag();
231}
232} // namespace
233
234double RooExpPoly::analyticalIntegral(int /*code*/, const char *rangeName) const
235{
236 const double xmin = _x.min(rangeName), xmax = _x.max(rangeName);
237 const unsigned sz = _coefList.size();
238 if (!sz)
239 return xmax - xmin;
240
241 std::vector<double> coefs;
242 coefs.reserve(sz);
243 const RooArgSet *nset = _coefList.nset();
244 for (auto c : _coefList) {
245 coefs.push_back(static_cast<RooAbsReal *>(c)->getVal(nset));
246 }
247
248 switch (_coefList.size() + _lowestOrder) {
249 case 1: return xmax - xmin;
250 case 2: {
251 const double a = coefs[1 - _lowestOrder];
252 if (a != 0) {
253 return 1. / a * (exp(a * xmax) - exp(a * xmin)) * (_lowestOrder == 0 ? exp(coefs[0]) : 1);
254 } else {
255 return xmax - xmin;
256 }
257 }
258 case 3: {
259 const double a = coefs[2 - _lowestOrder];
260 const double b = _lowestOrder == 2 ? 0. : coefs[1 - _lowestOrder];
261 const double c = _lowestOrder == 0 ? coefs[0] : 0.;
262 const double absa = std::abs(a);
263 const double sqrta = std::sqrt(absa);
264 if (a < 0) {
265 double d = ::deltaerf((-b + 2 * absa * xmax) / (2 * sqrta), (-b + 2 * absa * xmin) / (2 * sqrta));
266 double retval = exp(b * b / (4 * absa) + c) * std::sqrt(PI) * d / (2 * sqrta);
267 return retval;
268 } else if (a > 0) {
269 double d = ::deltaerfi((b + 2 * absa * xmax) / (2 * sqrta), (b + 2 * absa * xmin) / (2 * sqrta));
270 double retval = exp(-b * b / (4 * absa) + c) * std::sqrt(PI) * d / (2 * sqrta);
271 return retval;
272 } else if (b != 0) {
273 return 1. / b * (std::exp(b * xmax) - exp(b * xmin)) * exp(c);
274 } else {
275 return xmax - xmin;
276 }
277 }
278 }
279 return 0.;
280}
281
282////////////////////////////////////////////////////////////////////////////////
283
284std::string RooExpPoly::getFormulaExpression(bool expand) const
285{
286 std::stringstream ss;
287 ss << "exp(";
288 int order = _lowestOrder;
289 for (auto coef : _coefList) {
290 if (order != _lowestOrder)
291 ss << "+";
292 if (expand)
293 ss << ((RooAbsReal *)coef)->getVal();
294 else
295 ss << coef->GetName();
296 ss << "*pow(" << _x.GetName() << "," << order << ")";
297 ++order;
298 }
299 ss << ")";
300 return ss.str();
301}
#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:377
#define PI
#define R__ASSERT(e)
Definition TError.h:117
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
Storage_t::size_type size() const
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
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:62
double getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
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:55
bool add(const RooAbsArg &var, bool valueServer, bool shapeServer, bool silent)
Overloaded RooCollection_t::add() method insert object into set and registers object as server to own...
RooExpPoly implements a polynomial PDF of the form.
Definition RooExpPoly.h:18
double evaluate() const override
Evaluation.
double evaluateLog() const
std::string getFormulaExpression(bool expand) const
int lowestOrder() const
Return the order for the first coefficient in the list.
Definition RooExpPoly.h:33
double analyticalIntegral(int code, const char *rangeName=nullptr) const override
Implements the actual analytical integral(s) advertised by getAnalyticalIntegral.
RooListProxy _coefList
Definition RooExpPoly.h:46
int getAnalyticalIntegral(RooArgSet &allVars, RooArgSet &analVars, const char *rangeName=nullptr) const override
Interface function getAnalyticalIntergral advertises the analytical integrals that are supported.
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...
void adjustLimits()
RooRealProxy _x
Definition RooExpPoly.h:45
int _lowestOrder
Definition RooExpPoly.h:47
RooAbsReal const & x() const
Get the x variable.
Definition RooExpPoly.h:27
RooArgList const & coefList() const
Get the coefficient list.
Definition RooExpPoly.h:30
static std::complex< double > faddeeva(std::complex< double > z)
evaluate Faddeeva function for complex argument
Definition RooMath.cxx:30
RooRealVar represents a variable that can be changed from the outside.
Definition RooRealVar.h:40
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