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)
 
  139        std::string(
"Lambda needs to be negative when ") + 
_zeta.
GetName() + 
" is zero.");
 
  142#ifndef R__HAS_MATHMORE 
  143  throw std::logic_error(
"RooHypatia2 needs ROOT with mathmore enabled to access gsl functions.");
 
  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)
 
  165#ifndef R__HAS_MATHMORE 
  166  throw std::logic_error(
"RooHypatia2 needs ROOT with mathmore enabled to access gsl functions.");
 
  172const double logsq2pi = std::log(std::sqrt(
TMath::TwoPi()));
 
  173const double ln2 = std::log(2.);
 
  175double low_x_BK(
double nu, 
double x){
 
  176  return TMath::Gamma(nu)*std::pow(2., nu-1.)*std::pow(
x, -nu);
 
  179double low_x_LnBK(
double nu, 
double x){
 
  180  return std::log(
TMath::Gamma(nu)) + ln2*(nu-1.) - std::log(
x) * nu;
 
  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);
 
  190#ifdef R__HAS_MATHMORE 
  193  return std::numeric_limits<double>::signaling_NaN();
 
  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);
 
  205#ifdef R__HAS_MATHMORE 
  208  return std::numeric_limits<double>::signaling_NaN();
 
  213double LogEval(
double d, 
double l, 
double alpha, 
double beta, 
double delta) {
 
  214  const double gamma = alpha;
 
  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);
 
  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)) );
 
  226double diff_eval(
double d, 
double l, 
double alpha, 
double beta, 
double delta){
 
  227  const double gamma = alpha;
 
  228  const double dg = delta*
gamma;
 
  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;
 
  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;
 
  261  const double cons0 = std::sqrt(
_zeta);
 
  264  const double beta = 
_beta;
 
  271    const double cons1 = 
_sigma/std::sqrt(phi);
 
  272    const double alpha  = cons0/cons1;
 
  273    const double delta = cons0*cons1;
 
  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);
 
  281      out = A * std::pow(B-
d, -
_n);
 
  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);
 
  289      out = A * std::pow(B+
d, -
_n2);
 
  292      out = LogEval(
d, 
_lambda, alpha, beta, delta);
 
  295  else if (
_zeta < 0.) {
 
  296    coutE(Eval) << 
"The parameter " << 
_zeta.
GetName() << 
" of the RooHypatia2 " << 
GetName() << 
" cannot be < 0." << std::endl;
 
  299    const double delta = 
_sigma;
 
  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);
 
  309      out = A*std::pow(B-
d, -
_n);
 
  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);
 
  319      out =  A*std::pow(B+
d, -
_n2);
 
  322      out = std::exp(beta*
d) * std::pow(1. + 
d*
d/(delta*delta), 
_lambda - 0.5);
 
  326    coutE(Eval) << 
"zeta = 0 only supported for lambda < 0. lambda = " << 
double(
_lambda) << std::endl;
 
  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) {
 
  340  const bool zetaIsAlwaysLargerZero = !zeta.isBatch() && zeta[0] > 0.;
 
  341  const bool zetaIsAlwaysZero = !zeta.isBatch() && zeta[0] == 0.;
 
  343  for (std::size_t i = 0; i < 
N; ++i) {
 
  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];
 
  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;
 
  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]);
 
  363        output[i] = A * std::pow(B-
d, -
n[i]);
 
  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]);
 
  371        output[i] = A * std::pow(B+
d, -n2[i]);
 
  374        output[i] = LogEval(
d, lambda[i], alpha, beta[i], delta);
 
  378    if ((!zetaIsAlwaysLargerZero && zetaIsAlwaysZero) || (zeta[i] == 0. && lambda[i] < 0.)) {
 
  379      const double delta = 
sigma[i];
 
  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]);
 
  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]);
 
  399        output[i] =  A*std::pow(B+
d, -n2[i]);
 
  402        output[i] = std::exp(beta[i]*
d) * std::pow(1. + 
d*
d/(delta*delta), lambda[i] - 0.5);
 
  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);
 
  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);
 
  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) {
 
  440  const double cons0 = std::sqrt(zeta);
 
  441  const double asigma = 
a*
sigma;
 
  442  const double a2sigma = a2*
sigma;
 
  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;
 
  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);
 
  453    for (std::size_t i = 0; i < 
N; ++i) {
 
  454      const double d = 
x[i] - mu[i];
 
  457        output[i] = AB_l.first * std::pow(AB_l.second - 
d, -
n);
 
  459      else if (
d>a2sigma) {
 
  460        output[i] = AB_r.first * std::pow(AB_r.second + 
d, -n2);
 
  463        output[i] = LogEval(
d, lambda, alpha, beta, delta);
 
  466  } 
else if (zeta == 0. && lambda < 0.) {
 
  467    const double delta = 
sigma;
 
  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);
 
  472    for (std::size_t i = 0; i < 
N; ++i) {
 
  473      const double d = 
x[i] - mu[i];
 
  476        output[i] = AB_l.first*std::pow(AB_l.second - 
d, -
n);
 
  478      else if (
d > a2sigma) {
 
  479        output[i] = AB_r.first * std::pow(AB_r.second + 
d, -n2);
 
  482        output[i] = std::exp(beta*
d) * std::pow(1. + 
d*
d/(delta*delta), lambda - 0.5);
 
  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());
 
  512  if (
x.size()>1 && paramSizeSum==9) {
 
Bool_t isConstant() const
Check if the "Constant" attribute is set.
 
RooAbsReal is the common abstract base class for objects that represent a real value and implements f...
 
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.
 
RooArgSet is a container object that can hold multiple RooAbsArg objects.
 
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....
 
Double_t evaluate() const override
Evaluate this PDF / function / constant. Needs to be overridden by all derived classes.
 
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.
 
virtual const char * GetName() const
Returns name of object.
 
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...
 
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.
 
constexpr Double_t TwoPi()
 
This struct enables passing computation data around between elements of a computation graph.
 
RooSpan< double > makeBatch(const RooAbsReal *owner, std::size_t size)
Create a writable batch.
 
static void output(int code)