Logo ROOT  
Reference Guide
 
Loading...
Searching...
No Matches
RooHypatia2.cxx
Go to the documentation of this file.
1// Author: Stephan Hageboeck, CERN, Oct 2019
2// Based on RooIpatia2 by Diego Martinez Santos, Nikhef, Diego.Martinez.Santos@cern.ch
3/*****************************************************************************
4 * Project: RooFit *
5 * Authors: *
6 * WV, Wouter Verkerke, UC Santa Barbara, verkerke@slac.stanford.edu *
7 * DK, David Kirkby, UC Irvine, dkirkby@uci.edu *
8 * *
9 * Copyright (c) 2000-2019, Regents of the University of California *
10 * and Stanford University. All rights reserved. *
11 * *
12 * Redistribution and use in source and binary forms, *
13 * with or without modification, are permitted according to the terms *
14 * listed in LICENSE (http://roofit.sourceforge.net/license.txt) *
15 *****************************************************************************/
16
17/**
18 * \class RooHypatia2
19 *
20 * RooHypatia2 is the two-sided version of the Hypatia distribution described in https://arxiv.org/abs/1312.5000.
21 * \image html RooHypatia2.png
22 *
23 * It has a hyperbolic core of a crystal-ball-like function \f$ G \f$ and two tails:
24 * \f[
25 * \mathrm{Hypatia2}(x, \mu, \sigma, \lambda, \zeta, \beta, a_l, n_l, a_r, n_r) =
26 * \begin{cases}
27 * \frac{ G(\mu - a_l \sigma, \mu, \sigma, \lambda, \zeta, \beta) }
28 * { \left( 1 - \frac{x}{n_l G(\ldots)/G'(\ldots) - a_l\sigma } \right)^{n_l} }
29 * & \text{if } \frac{x-\mu}{\sigma} < -a_l \\
30 * \left( (x-\mu)^2 + A^2_\lambda(\zeta)\sigma^2 \right)^{\frac{1}{2}\lambda-\frac{1}{4}} e^{\beta(x-\mu)} K_{\lambda-\frac{1}{2}}
31 * \left( \zeta \sqrt{1+\left( \frac{x-\mu}{A_\lambda(\zeta)\sigma} \right)^2 } \right) \equiv G(x, \mu, \ldots)
32 * & \text{otherwise} \\
33 * \frac{ G(\mu + a_r \sigma, \mu, \sigma, \lambda, \zeta, \beta) }
34 * { \left( 1 - \frac{x}{-n_r G(\ldots)/G'(\ldots) - a_r\sigma } \right)^{n_r} }
35 * & \text{if } \frac{x-\mu}{\sigma} > a_r \\
36 * \end{cases}
37 * \f]
38 * Here, \f$ K_\lambda \f$ are the modified Bessel functions of the second kind
39 * ("irregular modified cylindrical Bessel functions" from the gsl,
40 * "special Bessel functions of the third kind"),
41 * and \f$ A^2_\lambda(\zeta) \f$ is a ratio of these:
42 * \f[
43 * A_\lambda^{2}(\zeta) = \frac{\zeta K_\lambda(\zeta)}{K_{\lambda+1}(\zeta)}
44 * \f]
45 *
46 * \if false
47 * TODO Enable once analytic integrals work.
48 * ### Analytical Integration
49 * The Hypatia distribution can be integrated analytically if \f$ \beta = \zeta = 0 \f$ and
50 * \f$ \lambda < 0 \f$. An analytic integral will only be used, though, if the parameters are **constant**
51 * at zero, and if \f$ \lambda < 0 \f$. This can be ensured as follows:
52 * ```
53 * RooRealVar beta("beta", "beta", 0.); // NOT beta("beta", "beta", 0., -1., 1.) This would allow it to float.
54 * RooRealVar zeta("zeta", "zeta", 0.);
55 * RooRealVar lambda("lambda", "lambda", -1., -10., -0.00001);
56 * ```
57 * In all other cases, slower / less accurate numeric integration will be used.
58 * Note that including `0.` in the value range of lambda excludes using analytic integrals.
59 * \endif
60 *
61 * ### Concavity
62 * Note that unless the parameters \f$ a_l,\ a_r \f$ are very large, the function has non-hyperbolic tails. This requires
63 * \f$ G \f$ to be strictly concave, *i.e.*, peaked, as otherwise the tails would yield imaginary numbers. Choosing \f$ \lambda,
64 * \beta, \zeta \f$ inappropriately will therefore lead to evaluation errors.
65 *
66 * Further, the original paper establishes that to keep the tails from rising,
67 * \f[
68 * \begin{split}
69 * \beta^2 &< \alpha^2 \\
70 * \Leftrightarrow \beta^2 &< \frac{\zeta^2}{\delta^2} = \frac{\zeta^2}{\sigma^2 A_{\lambda}^2(\zeta)}
71 * \end{split}
72 * \f]
73 * needs to be satisfied, unless the fit range is very restricted, because otherwise, the function rises in the tails.
74 *
75 *
76 * In case of evaluation errors, it is advisable to choose very large values for \f$ a_l,\ a_r \f$, tweak the parameters of the core region to
77 * make it concave, and re-enable the tails. Especially \f$ \beta \f$ needs to be close to zero.
78 *
79 * ## Relation to RooIpatia2
80 * This implementation is largely based on RooIpatia2, https://gitlab.cern.ch/lhcb/Urania/blob/master/PhysFit/P2VV/src/RooIpatia2.cxx,
81 * but there are differences:
82 * - At construction time, the Hypatia implementation checks if the range of parameters extends into regions where
83 * the function might be undefined.
84 * - Hypatia supports I/O to ROOT files.
85 * - Hypatia will support faster batched function evaluations.
86 * - Hypatia might support analytical integration for the case \f$ \zeta = \beta = 0, \lambda < 1 \f$.
87 *
88 * Because of these differences, and to avoid name clashes for setups where RooFit is used in an environment that also has
89 * RooIpatia2, class names need to be different.
90 */
91
92#include "RooHypatia2.h"
93#include "RooBatchCompute.h"
94#include "RooAbsReal.h"
95#include "RooHelpers.h"
96
97#include "TMath.h"
98#include "Math/SpecFunc.h"
99#include "ROOT/RConfig.hxx"
100
101#include <cmath>
102#include <algorithm>
103
104///////////////////////////////////////////////////////////////////////////////////////////
105/// Construct a new Hypatia2 PDF.
106/// \param[in] name Name of this instance.
107/// \param[in] title Title (for plotting)
108/// \param[in] x The variable of this distribution
109/// \param[in] lambda Shape parameter. Note that \f$ \lambda < 0 \f$ is required if \f$ \zeta = 0 \f$.
110/// \param[in] zeta Shape parameter (\f$ \zeta >= 0 \f$).
111/// \param[in] beta Asymmetry parameter \f$ \beta \f$. Symmetric case is \f$ \beta = 0 \f$,
112/// choose values close to zero.
113/// \param[in] sigma Width parameter. If \f$ \beta = 0, \ \sigma \f$ is the RMS width.
114/// \param[in] mu Location parameter. Shifts the distribution left/right.
115/// \param[in] a Start of the left tail (\f$ a \geq 0 \f$, to the left of the peak). Note that when setting \f$ a = \sigma = 1 \f$,
116/// the tail region is to the left of \f$ x = \mu - 1 \f$, so a should be positive.
117/// \param[in] n Shape parameter of left tail (\f$ n \ge 0 \f$). With \f$ n = 0 \f$, the function is constant.
118/// \param[in] a2 Start of right tail.
119/// \param[in] n2 Shape parameter of right tail (\f$ n2 \ge 0 \f$). With \f$ n2 = 0 \f$, the function is constant.
120RooHypatia2::RooHypatia2(const char *name, const char *title, RooAbsReal& x, RooAbsReal& lambda,
121 RooAbsReal& zeta, RooAbsReal& beta, RooAbsReal& argSigma, RooAbsReal& mu, RooAbsReal& a,
122 RooAbsReal& n, RooAbsReal& a2, RooAbsReal& n2) :
123 RooAbsPdf(name, title),
124 _x("x", "x", this, x),
125 _lambda("lambda", "Lambda", this, lambda),
126 _zeta("zeta", "zeta", this, zeta),
127 _beta("beta", "Asymmetry parameter beta", this, beta),
128 _sigma("sigma", "Width parameter sigma", this, argSigma),
129 _mu("mu", "Location parameter mu", this, mu),
130 _a("a", "Left tail location a", this, a),
131 _n("n", "Left tail parameter n", this, n),
132 _a2("a2", "Right tail location a2", this, a2),
133 _n2("n2", "Right tail parameter n2", this, n2)
134{
135 RooHelpers::checkRangeOfParameters(this, {&argSigma}, 0.);
136 RooHelpers::checkRangeOfParameters(this, {&zeta, &n, &n2, &a, &a2}, 0., std::numeric_limits<double>::max(), true);
137 if (zeta.getVal() == 0. && zeta.isConstant()) {
138 RooHelpers::checkRangeOfParameters(this, {&lambda}, -std::numeric_limits<double>::max(), 0., false,
139 std::string("Lambda needs to be negative when ") + _zeta.GetName() + " is zero.");
140 }
141
142#ifndef R__HAS_MATHMORE
143 throw std::logic_error("RooHypatia2 needs ROOT with mathmore enabled to access gsl functions.");
144#endif
145}
146
147
148///////////////////////////////////////////////////////////////////////////////////////////
149/// Copy a new Hypatia2 PDF.
150/// \param[in] other Original to copy from.
151/// \param[in] name Optional new name.
152RooHypatia2::RooHypatia2(const RooHypatia2& other, const char* name) :
153 RooAbsPdf(other, name),
154 _x("x", this, other._x),
155 _lambda("lambda", this, other._lambda),
156 _zeta("zeta", this, other._zeta),
157 _beta("beta", this, other._beta),
158 _sigma("sigma", this, other._sigma),
159 _mu("mu", this, other._mu),
160 _a("a", this, other._a),
161 _n("n", this, other._n),
162 _a2("a2", this, other._a2),
163 _n2("n2", this, other._n2)
164{
165#ifndef R__HAS_MATHMORE
166 throw std::logic_error("RooHypatia2 needs ROOT with mathmore enabled to access gsl functions.");
167#endif
168}
169
170namespace {
171const double sq2pi_inv = 1./std::sqrt(TMath::TwoPi());
172const double logsq2pi = std::log(std::sqrt(TMath::TwoPi()));
173const double ln2 = std::log(2.);
174
175double low_x_BK(double nu, double x){
176 return TMath::Gamma(nu)*std::pow(2., nu-1.)*std::pow(x, -nu);
177}
178
179double low_x_LnBK(double nu, double x){
180 return std::log(TMath::Gamma(nu)) + ln2*(nu-1.) - std::log(x) * nu;
181}
182
183double besselK(double ni, double x) {
184 const double nu = std::fabs(ni);
185 if ((x < 1.e-06 && nu > 0.) ||
186 (x < 1.e-04 && nu > 0. && nu < 55.) ||
187 (x < 0.1 && nu >= 55.) )
188 return low_x_BK(nu, x);
189
190#ifdef R__HAS_MATHMORE
191 return ROOT::Math::cyl_bessel_k(nu, x);
192#else
193 return std::numeric_limits<double>::signaling_NaN();
194#endif
195
196}
197
198double LnBesselK(double ni, double x) {
199 const double nu = std::fabs(ni);
200 if ((x < 1.e-06 && nu > 0.) ||
201 (x < 1.e-04 && nu > 0. && nu < 55.) ||
202 (x < 0.1 && nu >= 55.) )
203 return low_x_LnBK(nu, x);
204
205#ifdef R__HAS_MATHMORE
206 return std::log(ROOT::Math::cyl_bessel_k(nu, x));
207#else
208 return std::numeric_limits<double>::signaling_NaN();
209#endif
210}
211
212
213double LogEval(double d, double l, double alpha, double beta, double delta) {
214 const double gamma = alpha;//std::sqrt(alpha*alpha-beta*beta);
215 const double dg = delta*gamma;
216 const double thing = delta*delta + d*d;
217 const double logno = l*std::log(gamma/delta) - logsq2pi - LnBesselK(l, dg);
218
219 return std::exp(logno + beta*d
220 + (0.5-l)*(std::log(alpha)-0.5*std::log(thing))
221 + LnBesselK(l-0.5, alpha*std::sqrt(thing)) );// + std::log(std::fabs(beta)+0.0001) );
222
223}
224
225
226double diff_eval(double d, double l, double alpha, double beta, double delta){
227 const double gamma = alpha;
228 const double dg = delta*gamma;
229
230 const double thing = delta*delta + d*d;
231 const double sqrthing = std::sqrt(thing);
232 const double alphasq = alpha*sqrthing;
233 const double no = std::pow(gamma/delta, l)/besselK(l, dg)*sq2pi_inv;
234 const double ns1 = 0.5-l;
235
236 return no * std::pow(alpha, ns1) * std::pow(thing, l/2. - 1.25)
237 * ( -d * alphasq * (besselK(l - 1.5, alphasq)
238 + besselK(l + 0.5, alphasq))
239 + (2.*(beta*thing + d*l) - d) * besselK(ns1, alphasq) )
240 * std::exp(beta*d) * 0.5;
241}
242
243/*
244double Gauss2F1(double a, double b, double c, double x){
245 if (fabs(x) <= 1.) {
246 return ROOT::Math::hyperg(a, b, c, x);
247 } else {
248 return ROOT::Math::hyperg(c-a, b, c, 1-1/(1-x))/std::pow(1-x, b);
249 }
250}
251
252double stIntegral(double d1, double delta, double l){
253 return d1 * Gauss2F1(0.5, 0.5-l, 3./2, -d1*d1/(delta*delta));
254}
255*/
256}
257
259{
260 const double d = _x-_mu;
261 const double cons0 = std::sqrt(_zeta);
262 const double asigma = _a*_sigma;
263 const double a2sigma = _a2*_sigma;
264 const double beta = _beta;
265 double out = 0.;
266
267 if (_zeta > 0.) {
268 // careful if zeta -> 0. You can implement a function for the ratio,
269 // but careful again that |nu + 1 | != |nu| + 1 so you have to deal with the signs
270 const double phi = besselK(_lambda+1., _zeta) / besselK(_lambda, _zeta);
271 const double cons1 = _sigma/std::sqrt(phi);
272 const double alpha = cons0/cons1;
273 const double delta = cons0*cons1;
274
275 if (d < -asigma){
276 const double k1 = LogEval(-asigma, _lambda, alpha, beta, delta);
277 const double k2 = diff_eval(-asigma, _lambda, alpha, beta, delta);
278 const double B = -asigma + _n*k1/k2;
279 const double A = k1 * std::pow(B+asigma, _n);
280
281 out = A * std::pow(B-d, -_n);
282 }
283 else if (d>a2sigma) {
284 const double k1 = LogEval(a2sigma, _lambda, alpha, beta, delta);
285 const double k2 = diff_eval(a2sigma, _lambda, alpha, beta, delta);
286 const double B = -a2sigma - _n2*k1/k2;
287 const double A = k1 * std::pow(B+a2sigma, _n2);
288
289 out = A * std::pow(B+d, -_n2);
290 }
291 else {
292 out = LogEval(d, _lambda, alpha, beta, delta);
293 }
294 }
295 else if (_zeta < 0.) {
296 coutE(Eval) << "The parameter " << _zeta.GetName() << " of the RooHypatia2 " << GetName() << " cannot be < 0." << std::endl;
297 }
298 else if (_lambda < 0.) {
299 const double delta = _sigma;
300
301 if (d < -asigma ) {
302 const double cons1 = std::exp(-beta*asigma);
303 const double phi = 1. + _a * _a;
304 const double k1 = cons1 * std::pow(phi, _lambda-0.5);
305 const double k2 = beta*k1 - cons1*(_lambda-0.5) * std::pow(phi, _lambda-1.5) * 2.*_a/delta;
306 const double B = -asigma + _n*k1/k2;
307 const double A = k1*std::pow(B+asigma, _n);
308
309 out = A*std::pow(B-d, -_n);
310 }
311 else if (d > a2sigma) {
312 const double cons1 = std::exp(beta*a2sigma);
313 const double phi = 1. + _a2*_a2;
314 const double k1 = cons1 * std::pow(phi, _lambda-0.5);
315 const double k2 = beta*k1 + cons1*(_lambda-0.5) * std::pow(phi, _lambda-1.5) * 2.*_a2/delta;
316 const double B = -a2sigma - _n2*k1/k2;
317 const double A = k1*std::pow(B+a2sigma, _n2);
318
319 out = A*std::pow(B+d, -_n2);
320 }
321 else {
322 out = std::exp(beta*d) * std::pow(1. + d*d/(delta*delta), _lambda - 0.5);
323 }
324 }
325 else {
326 coutE(Eval) << "zeta = 0 only supported for lambda < 0. lambda = " << double(_lambda) << std::endl;
327 }
328
329 return out;
330}
331
332namespace {
333//////////////////////////////////////////////////////////////////////////////////////////
334/// The generic compute function that recalculates everything for every loop iteration.
335/// This is only needed in the rare case where a parameter is used as an observable.
336template<typename Tx, typename Tl, typename Tz, typename Tb, typename Ts, typename Tm, typename Ta, typename Tn,
337typename Ta2, typename Tn2>
338void compute(RooSpan<double> output, Tx x, Tl lambda, Tz zeta, Tb beta, Ts sigma, Tm mu, Ta a, Tn n, Ta2 a2, Tn2 n2) {
339 const auto N = output.size();
340 const bool zetaIsAlwaysLargerZero = !zeta.isBatch() && zeta[0] > 0.;
341 const bool zetaIsAlwaysZero = !zeta.isBatch() && zeta[0] == 0.;
342
343 for (std::size_t i = 0; i < N; ++i) {
344
345 const double d = x[i] - mu[i];
346 const double cons0 = std::sqrt(zeta[i]);
347 const double asigma = a[i]*sigma[i];
348 const double a2sigma = a2[i]*sigma[i];
349// const double beta = beta[i];
350
351 if (zetaIsAlwaysLargerZero || zeta[i] > 0.) {
352 const double phi = besselK(lambda[i]+1., zeta[i]) / besselK(lambda[i], zeta[i]);
353 const double cons1 = sigma[i]/std::sqrt(phi);
354 const double alpha = cons0/cons1;
355 const double delta = cons0*cons1;
356
357 if (d < -asigma){
358 const double k1 = LogEval(-asigma, lambda[i], alpha, beta[i], delta);
359 const double k2 = diff_eval(-asigma, lambda[i], alpha, beta[i], delta);
360 const double B = -asigma + n[i]*k1/k2;
361 const double A = k1 * std::pow(B+asigma, n[i]);
362
363 output[i] = A * std::pow(B-d, -n[i]);
364 }
365 else if (d>a2sigma) {
366 const double k1 = LogEval(a2sigma, lambda[i], alpha, beta[i], delta);
367 const double k2 = diff_eval(a2sigma, lambda[i], alpha, beta[i], delta);
368 const double B = -a2sigma - n2[i]*k1/k2;
369 const double A = k1 * std::pow(B+a2sigma, n2[i]);
370
371 output[i] = A * std::pow(B+d, -n2[i]);
372 }
373 else {
374 output[i] = LogEval(d, lambda[i], alpha, beta[i], delta);
375 }
376 }
377
378 if ((!zetaIsAlwaysLargerZero && zetaIsAlwaysZero) || (zeta[i] == 0. && lambda[i] < 0.)) {
379 const double delta = sigma[i];
380
381 if (d < -asigma ) {
382 const double cons1 = std::exp(-beta[i]*asigma);
383 const double phi = 1. + a[i] * a[i];
384 const double k1 = cons1 * std::pow(phi, lambda[i]-0.5);
385 const double k2 = beta[i]*k1 - cons1*(lambda[i]-0.5) * std::pow(phi, lambda[i]-1.5) * 2.*a[i]/delta;
386 const double B = -asigma + n[i]*k1/k2;
387 const double A = k1*std::pow(B+asigma, n[i]);
388
389 output[i] = A*std::pow(B-d, -n[i]);
390 }
391 else if (d > a2sigma) {
392 const double cons1 = std::exp(beta[i]*a2sigma);
393 const double phi = 1. + a2[i]*a2[i];
394 const double k1 = cons1 * std::pow(phi, lambda[i]-0.5);
395 const double k2 = beta[i]*k1 + cons1*(lambda[i]-0.5) * std::pow(phi, lambda[i]-1.5) * 2.*a2[i]/delta;
396 const double B = -a2sigma - n2[i]*k1/k2;
397 const double A = k1*std::pow(B+a2sigma, n2[i]);
398
399 output[i] = A*std::pow(B+d, -n2[i]);
400 }
401 else {
402 output[i] = std::exp(beta[i]*d) * std::pow(1. + d*d/(delta*delta), lambda[i] - 0.5);
403 }
404 }
405 }
406}
407
408template<bool right>
409std::pair<double, double> computeAB_zetaNZero(double asigma, double lambda, double alpha, double beta, double delta, double n) {
410 const double k1 = LogEval(right ? asigma : -asigma, lambda, alpha, beta, delta);
411 const double k2 = diff_eval(right ? asigma : -asigma, lambda, alpha, beta, delta);
412 const double B = -asigma + (right ? -1 : 1.) * n*k1/k2;
413 const double A = k1 * std::pow(B+asigma, n);
414
415 return {A, B};
416}
417
418template<bool right>
419std::pair<double, double> computeAB_zetaZero(double beta, double asigma, double a, double lambda, double delta, double n) {
420 const double cons1 = std::exp((right ? 1. : -1.) * beta*asigma);
421 const double phi = 1. + a * a;
422 const double k1 = cons1 * std::pow(phi, lambda-0.5);
423 const double k2 = beta*k1 + (right ? 1. : -1) * cons1*(lambda-0.5) * std::pow(phi, lambda-1.5) * 2.*a/delta;
424 const double B = -asigma + (right ? -1. : 1.) * n*k1/k2;
425 const double A = k1*std::pow(B+asigma, n);
426
427 return {A, B};
428}
429
431//////////////////////////////////////////////////////////////////////////////////////////
432/// A specialised compute function where x is an observable, and all parameters are used as
433/// parameters. Since many things can be calculated outside of the loop, it is faster.
435 BracketAdapter<double> lambda, BracketAdapter<double> zeta, BracketAdapter<double> beta,
436 BracketAdapter<double> sigma, BracketAdapter<double> mu,
437 BracketAdapter<double> a, BracketAdapter<double> n, BracketAdapter<double> a2, BracketAdapter<double> n2) {
438 const auto N = output.size();
439
440 const double cons0 = std::sqrt(zeta);
441 const double asigma = a*sigma;
442 const double a2sigma = a2*sigma;
443
444 if (zeta > 0.) {
445 const double phi = besselK(lambda+1., zeta) / besselK(lambda, zeta);
446 const double cons1 = sigma/std::sqrt(phi);
447 const double alpha = cons0/cons1;
448 const double delta = cons0*cons1;
449
450 const auto AB_l = computeAB_zetaNZero<false>(asigma, lambda, alpha, beta, delta, n);
451 const auto AB_r = computeAB_zetaNZero<true>(a2sigma, lambda, alpha, beta, delta, n2);
452
453 for (std::size_t i = 0; i < N; ++i) {
454 const double d = x[i] - mu[i];
455
456 if (d < -asigma){
457 output[i] = AB_l.first * std::pow(AB_l.second - d, -n);
458 }
459 else if (d>a2sigma) {
460 output[i] = AB_r.first * std::pow(AB_r.second + d, -n2);
461 }
462 else {
463 output[i] = LogEval(d, lambda, alpha, beta, delta);
464 }
465 }
466 } else if (zeta == 0. && lambda < 0.) {
467 const double delta = sigma;
468
469 const auto AB_l = computeAB_zetaZero<false>(beta, asigma, a, lambda, delta, n);
470 const auto AB_r = computeAB_zetaZero<true>(beta, a2sigma, a2, lambda, delta, n2);
471
472 for (std::size_t i = 0; i < N; ++i) {
473 const double d = x[i] - mu[i];
474
475 if (d < -asigma ) {
476 output[i] = AB_l.first*std::pow(AB_l.second - d, -n);
477 }
478 else if (d > a2sigma) {
479 output[i] = AB_r.first * std::pow(AB_r.second + d, -n2);
480 }
481 else {
482 output[i] = std::exp(beta*d) * std::pow(1. + d*d/(delta*delta), lambda - 0.5);
483 }
484 }
485 }
486}
487
488}
489
491 using namespace RooBatchCompute;
492
493 auto x = _x->getValues(evalData, normSet);
494 auto lambda = _lambda->getValues(evalData, normSet);
495 auto zeta = _zeta->getValues(evalData, normSet);
496 auto beta = _beta->getValues(evalData, normSet);
497 auto sig = _sigma->getValues(evalData, normSet);
498 auto mu = _mu->getValues(evalData, normSet);
499 auto a = _a->getValues(evalData, normSet);
500 auto n = _n->getValues(evalData, normSet);
501 auto a2 = _a2->getValues(evalData, normSet);
502 auto n2 = _n2->getValues(evalData, normSet);
503
504 size_t paramSizeSum=0, batchSize = x.size();
505 for (const auto& i:{lambda, zeta, beta, sig, mu, a, n, a2, n2}) {
506 paramSizeSum += i.size();
507 batchSize = std::max(batchSize, i.size());
508 }
509 RooSpan<double> output = evalData.makeBatch(this, batchSize);
510
511 // Run high performance compute if only x has multiple values
512 if (x.size()>1 && paramSizeSum==9) {
513 compute(output, x,
518 } else {
525 }
526 return output;
527}
528
529
530/* Analytical integrals need testing.
531
532Int_t RooHypatia2::getAnalyticalIntegral(RooArgSet& allVars, RooArgSet& analVars, const char*) const
533{
534 if (matchArgs(allVars, analVars, _x)
535 && _beta == 0. && _beta.arg().isConstant()
536 && _zeta == 0. && _zeta.arg().isConstant()
537 && _lambda.max() < 0.) return 1;
538 return 0 ;
539}
540
541
542
543double RooHypatia2::analyticalIntegral(Int_t code, const char* rangeName) const
544{
545 if (_beta != 0. || _zeta != 0. || _lambda >= 0) {
546 auto& logstream = coutF(Integration) << "When the PDF " << GetName()
547 << " was constructed, beta,zeta were 0, lambda<0 and all three were constant.\n"
548 << "This allowed for analytic integration, but now, numeric integration would be required.\n"
549 << "These parameters must either be kept constant, or be floating from the beginning. "
550 << " Terminating fit ..."
551 << std::endl;
552 RooArgSet vars;
553 vars.add(_beta.arg());
554 vars.add(_zeta.arg());
555 vars.add(_lambda.arg());
556 vars.printStream(logstream, vars.defaultPrintContents(nullptr), RooPrintable::kVerbose);
557 throw std::runtime_error("RooHypatia2 cannot be integrated analytically since parameters changed.");
558 }
559
560 // The formulae to follow still use beta and zeta to facilitate comparisons with the
561 // evaluate function. Since beta and zeta are zero, all relevant terms will be disabled
562 // by defining these two constexpr:
563 constexpr double beta = 0.;
564 constexpr double cons1 = 1.;
565
566 if (code != 1) {
567 throw std::logic_error("Trying to compute analytic integral with unknown configuration.");
568 }
569
570 const double asigma = _a * _sigma;
571 const double a2sigma = _a2 * _sigma;
572 const double d0 = _x.min(rangeName) - _mu;
573 const double d1 = _x.max(rangeName) - _mu;
574
575
576 double delta;
577 if (_lambda <= -1.0) {
578 delta = _sigma * sqrt(-2. + -2.*_lambda);
579 }
580 else {
581 delta = _sigma;
582 }
583 const double deltaSq = delta*delta;
584
585 if ((d0 > -asigma) && (d1 < a2sigma)) {
586 return stIntegral(d1, delta, _lambda) - stIntegral(d0, delta, _lambda);
587 }
588
589 if (d0 > a2sigma) {
590 const double phi = 1. + a2sigma*a2sigma/deltaSq;
591 const double k1 = cons1*std::pow(phi,_lambda-0.5);
592 const double k2 = beta*k1+ cons1*(_lambda-0.5)*std::pow(phi,_lambda-1.5)*2.*a2sigma/deltaSq;
593 const double B = -a2sigma - _n2*k1/k2;
594 const double A = k1*std::pow(B+a2sigma,_n2);
595 return A*(std::pow(B+d1,1-_n2)/(1-_n2) -std::pow(B+d0,1-_n2)/(1-_n2) ) ;
596
597 }
598
599 if (d1 < -asigma) {
600 const double phi = 1. + asigma*asigma/deltaSq;
601 const double k1 = cons1*std::pow(phi,_lambda-0.5);
602 const double k2 = beta*k1- cons1*(_lambda-0.5)*std::pow(phi,_lambda-1.5)*2*asigma/deltaSq;
603 const double B = -asigma + _n*k1/k2;
604 const double A = k1*std::pow(B+asigma,_n);
605 const double I0 = A*std::pow(B-d0,1-_n)/(_n-1);
606 const double I1 = A*std::pow(B-d1,1-_n)/(_n-1);
607 return I1 - I0;
608 }
609
610
611 double I0;
612 double I1a = 0;
613 double I1b = 0;
614 if (d0 <-asigma) {
615 const double phi = 1. + asigma*asigma/deltaSq;
616 const double k1 = cons1*std::pow(phi,_lambda-0.5);
617 const double k2 = beta*k1- cons1*(_lambda-0.5)*std::pow(phi,_lambda-1.5)*2*asigma/deltaSq;
618 const double B = -asigma + _n*k1/k2;
619 const double A = k1*std::pow(B+asigma,_n);
620 I0 = A*std::pow(B-d0,1-_n)/(_n-1);
621 I1a = A*std::pow(B+asigma,1-_n)/(_n-1) - stIntegral(-asigma, delta, _lambda);
622 }
623 else {
624 I0 = stIntegral(d0, delta, _lambda);
625 }
626
627 if (d1 > a2sigma) {
628 const double phi = 1. + a2sigma*a2sigma/deltaSq;
629 const double k1 = cons1*std::pow(phi,_lambda-0.5);
630 const double k2 = beta*k1+ cons1*(_lambda-0.5)*std::pow(phi,_lambda-1.5)*2.*a2sigma/deltaSq;
631 const double B = -a2sigma - _n2*k1/k2;
632 const double A = k1*std::pow(B+a2sigma,_n2);
633 I1b = A*(std::pow(B+d1,1-_n2)/(1-_n2) -std::pow(B+a2sigma,1-_n2)/(1-_n2) ) - stIntegral(d1, delta, _lambda) + stIntegral(a2sigma,delta, _lambda) ;
634 }
635
636 const double I1 = stIntegral(d1, delta, _lambda) + I1a + I1b;
637 return I1 - I0;
638}
639
640*/
double
#define d(i)
Definition RSha256.hxx:102
#define a(i)
Definition RSha256.hxx:99
#define coutE(a)
#define N
char name[80]
Definition TGX11.cxx:110
Bool_t isConstant() const
Check if the "Constant" attribute is set.
Definition RooAbsArg.h:380
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
Definition RooAbsReal.h:61
virtual RooSpan< const double > getValues(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet=nullptr) const
by this change, please consult the release notes for ROOT 6.24 for guidance on how to make this trans...
Double_t getVal(const RooArgSet *normalisationSet=nullptr) const
Evaluate object.
Definition RooAbsReal.h:91
RooArgSet is a container object that can hold multiple RooAbsArg objects.
Definition RooArgSet.h:29
Little adapter that gives a bracket operator to types that don't have one.
RooHypatia2 is the two-sided version of the Hypatia distribution described in https://arxiv....
Definition RooHypatia2.h:25
RooRealProxy _beta
Definition RooHypatia2.h:47
RooRealProxy _n
Definition RooHypatia2.h:51
RooRealProxy _a
Definition RooHypatia2.h:50
RooRealProxy _zeta
Definition RooHypatia2.h:46
RooRealProxy _n2
Definition RooHypatia2.h:53
RooRealProxy _a2
Definition RooHypatia2.h:52
RooRealProxy _mu
Definition RooHypatia2.h:49
RooRealProxy _sigma
Definition RooHypatia2.h:48
RooRealProxy _lambda
Definition RooHypatia2.h:45
Double_t evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
RooRealProxy _x
Definition RooHypatia2.h:44
RooSpan< double > evaluateSpan(RooBatchCompute::RunContext &evalData, const RooArgSet *normSet) const override
Evaluate this object for a batch/span of data points.
A simple container to hold a batch of data values.
Definition RooSpan.h:34
virtual const char * GetName() const
Returns name of object.
Definition TNamed.h:47
double beta(double x, double y)
Calculates the beta function.
double cyl_bessel_k(double nu, double x)
Calculates the modified Bessel functions of the second kind (also called irregular modified (cylindri...
const Double_t sigma
Double_t x[n]
Definition legend1.C:17
const Int_t n
Definition legend1.C:16
double gamma(double x)
Namespace for dispatching RooFit computations to various backends.
void checkRangeOfParameters(const RooAbsReal *callingClass, std::initializer_list< const RooAbsReal * > pars, double min=-std::numeric_limits< double >::max(), double max=std::numeric_limits< double >::max(), bool limitsInAllowedRange=false, std::string extraMessage="")
Check if the parameters have a range, and warn if the range extends below / above the set limits.
Double_t Gamma(Double_t z)
Computation of gamma(z) for all z.
Definition TMath.cxx:348
constexpr Double_t TwoPi()
Definition TMath.h:44
This struct enables passing computation data around between elements of a computation graph.
Definition RunContext.h:31
RooSpan< double > makeBatch(const RooAbsReal *owner, std::size_t size)
Create a writable batch.
auto * l
Definition textangle.C:4
static void output(int code)
Definition gifencode.c:226