File indexing completed on 2026-06-17 08:33:51
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014 #ifndef RooFit_Detail_MathFuncs_h
0015 #define RooFit_Detail_MathFuncs_h
0016
0017 #include <ROOT/RConfig.hxx> // for R__HAS_CLAD
0018
0019 #include <TMath.h>
0020 #include <Math/PdfFuncMathCore.h>
0021 #include <Math/ProbFuncMathCore.h>
0022
0023 #include <algorithm>
0024 #include <cmath>
0025 #include <stdexcept>
0026
0027 namespace RooFit::Detail::MathFuncs {
0028
0029
0030
0031 inline double binomial(int n, int k)
0032 {
0033 if (n < 0 || k < 0 || n < k)
0034 return TMath::SignalingNaN();
0035 if (k == 0 || n == k)
0036 return 1;
0037
0038 int k1 = std::min(k, n - k);
0039 int k2 = n - k1;
0040 double fact = k2 + 1;
0041 for (double i = k1; i > 1.; --i) {
0042 fact *= (k2 + i) / i;
0043 }
0044 return fact;
0045 }
0046
0047
0048 template <typename DoubleArray>
0049 double bernstein(double x, double xmin, double xmax, DoubleArray coefs, int nCoefs)
0050 {
0051 double xScaled = (x - xmin) / (xmax - xmin);
0052 int degree = nCoefs - 1;
0053
0054
0055 if (degree < 0) {
0056 return TMath::SignalingNaN();
0057 } else if (degree == 0) {
0058 return coefs[0];
0059 } else if (degree == 1) {
0060
0061 double a0 = coefs[0];
0062 double a1 = coefs[1] - a0;
0063 return a1 * xScaled + a0;
0064
0065 } else if (degree == 2) {
0066
0067 double a0 = coefs[0];
0068 double a1 = 2 * (coefs[1] - a0);
0069 double a2 = coefs[2] - a1 - a0;
0070 return (a2 * xScaled + a1) * xScaled + a0;
0071 }
0072
0073 double t = xScaled;
0074 double s = 1. - xScaled;
0075
0076 double result = coefs[0] * s;
0077 for (int i = 1; i < degree; i++) {
0078 result = (result + t * binomial(degree, i) * coefs[i]) * s;
0079 t *= xScaled;
0080 }
0081 result += t * coefs[degree];
0082
0083 return result;
0084 }
0085
0086
0087 inline double gaussian(double x, double mean, double sigma)
0088 {
0089 const double arg = x - mean;
0090 const double sig = sigma;
0091 return std::exp(-0.5 * arg * arg / (sig * sig));
0092 }
0093
0094 template <typename DoubleArray>
0095 double product(DoubleArray factors, std::size_t nFactors)
0096 {
0097 double out = 1.0;
0098 for (std::size_t i = 0; i < nFactors; ++i) {
0099 out *= factors[i];
0100 }
0101 return out;
0102 }
0103
0104
0105 inline double ratio(double numerator, double denominator)
0106 {
0107 return numerator / denominator;
0108 }
0109
0110 inline double bifurGauss(double x, double mean, double sigmaL, double sigmaR)
0111 {
0112
0113
0114 if (x < mean)
0115 return gaussian(x, mean, sigmaL);
0116 return gaussian(x, mean, sigmaR);
0117 }
0118
0119 inline double efficiency(double effFuncVal, int catIndex, int sigCatIndex)
0120 {
0121
0122 effFuncVal = std::clamp(effFuncVal, 0.0, 1.0);
0123
0124 if (catIndex == sigCatIndex)
0125 return effFuncVal;
0126 else
0127 return 1 - effFuncVal;
0128 }
0129
0130
0131 template <bool pdfMode = false, typename DoubleArray>
0132 double polynomial(DoubleArray coeffs, int nCoeffs, int lowestOrder, double x)
0133 {
0134 double retVal = coeffs[nCoeffs - 1];
0135 for (int i = nCoeffs - 2; i >= 0; i--) {
0136 retVal = coeffs[i] + x * retVal;
0137 }
0138 retVal = retVal * std::pow(x, lowestOrder);
0139 return retVal + (pdfMode && lowestOrder > 0 ? 1.0 : 0.0);
0140 }
0141
0142 template <typename DoubleArray>
0143 double chebychev(DoubleArray coeffs, unsigned int nCoeffs, double x_in, double xMin, double xMax)
0144 {
0145
0146 const double xPrime = (x_in - 0.5 * (xMax + xMin)) / (0.5 * (xMax - xMin));
0147
0148
0149 double sum = 1.;
0150 if (nCoeffs > 0) {
0151 double curr = xPrime;
0152 double twox = 2 * xPrime;
0153 double last = 1;
0154 double newval = twox * curr - last;
0155 last = curr;
0156 curr = newval;
0157 for (unsigned int i = 0; nCoeffs != i; ++i) {
0158 sum += last * coeffs[i];
0159 newval = twox * curr - last;
0160 last = curr;
0161 curr = newval;
0162 }
0163 }
0164 return sum;
0165 }
0166
0167 template <typename DoubleArray>
0168 double multipdf(int idx, DoubleArray pdfs)
0169 {
0170
0171
0172
0173
0174
0175 return pdfs[idx];
0176 }
0177
0178 template <typename DoubleArray>
0179 double constraintSum(DoubleArray comp, unsigned int compSize)
0180 {
0181 double sum = 0;
0182 #if defined(__CLING__) && defined(R__HAS_CLAD)
0183 #pragma clad checkpoint loop
0184 #endif
0185 for (unsigned int i = 0; i < compSize; i++) {
0186 sum -= std::log(comp[i]);
0187 }
0188 return sum;
0189 }
0190
0191 inline unsigned int uniformBinNumber(double low, double high, double val, unsigned int numBins, double coef)
0192 {
0193 double binWidth = (high - low) / numBins;
0194 return coef * (val >= high ? numBins - 1 : std::abs((val - low) / binWidth));
0195 }
0196
0197 template <typename DoubleArray>
0198 unsigned int rawBinNumber(double x, DoubleArray boundaries, std::size_t nBoundaries)
0199 {
0200 DoubleArray end = boundaries + nBoundaries;
0201 DoubleArray it = std::lower_bound(boundaries, end, x);
0202
0203 while (boundaries != it && (end == it || end == it + 1 || x < *it)) {
0204 --it;
0205 }
0206 return it - boundaries;
0207 }
0208
0209 template <typename DoubleArray>
0210 unsigned int binNumber(double x, double coef, DoubleArray boundaries, unsigned int nBoundaries, int nbins, int blo)
0211 {
0212 const int rawBin = rawBinNumber(x, boundaries, nBoundaries);
0213 int tmp = std::min(nbins, rawBin - blo);
0214 return coef * std::max(0, tmp);
0215 }
0216
0217 template <typename DoubleArray>
0218 double interpolate1d(double low, double high, double val, unsigned int numBins, DoubleArray vals)
0219 {
0220 double binWidth = (high - low) / numBins;
0221 int idx = val >= high ? numBins - 1 : std::abs((val - low) / binWidth);
0222
0223
0224 double central = low + (idx + 0.5) * binWidth;
0225 if (val > low + 0.5 * binWidth && val < high - 0.5 * binWidth) {
0226 double slope;
0227 if (val < central) {
0228 slope = vals[idx] - vals[idx - 1];
0229 } else {
0230 slope = vals[idx + 1] - vals[idx];
0231 }
0232 return vals[idx] + slope * (val - central) / binWidth;
0233 }
0234
0235 return vals[idx];
0236 }
0237
0238 inline double poisson(double x, double par)
0239 {
0240 if (par < 0)
0241 return TMath::QuietNaN();
0242
0243 if (x < 0) {
0244 return 0;
0245 } else if (x == 0.0) {
0246 return std::exp(-par);
0247 } else {
0248 double out = x * std::log(par) - TMath::LnGamma(x + 1.) - par;
0249 return std::exp(out);
0250 }
0251 }
0252
0253 inline double flexibleInterpSingle(unsigned int code, double low, double high, double boundary, double nominal,
0254 double paramVal, double res)
0255 {
0256 if (code == 0) {
0257
0258 if (paramVal > 0) {
0259 return paramVal * (high - nominal);
0260 } else {
0261 return paramVal * (nominal - low);
0262 }
0263 } else if (code == 1) {
0264
0265 if (paramVal >= 0) {
0266 return res * (std::pow(high / nominal, +paramVal) - 1);
0267 } else {
0268 return res * (std::pow(low / nominal, -paramVal) - 1);
0269 }
0270 } else if (code == 2) {
0271
0272 double a = 0.5 * (high + low) - nominal;
0273 double b = 0.5 * (high - low);
0274 double c = 0;
0275 if (paramVal > 1) {
0276 return (2 * a + b) * (paramVal - 1) + high - nominal;
0277 } else if (paramVal < -1) {
0278 return -1 * (2 * a - b) * (paramVal + 1) + low - nominal;
0279 } else {
0280 return a * paramVal * paramVal + b * paramVal + c;
0281 }
0282
0283
0284
0285
0286 } else if (code == 4 || code == 6) {
0287 double x = paramVal;
0288 double mod = 1.0;
0289 if (code == 6) {
0290 high /= nominal;
0291 low /= nominal;
0292 nominal = 1;
0293 }
0294 if (x >= boundary) {
0295 mod = x * (high - nominal);
0296 } else if (x <= -boundary) {
0297 mod = x * (nominal - low);
0298 } else {
0299
0300 double t = x / boundary;
0301 double eps_plus = high - nominal;
0302 double eps_minus = nominal - low;
0303 double S = 0.5 * (eps_plus + eps_minus);
0304 double A = 0.0625 * (eps_plus - eps_minus);
0305
0306 mod = x * (S + t * A * (15 + t * t * (-10 + t * t * 3)));
0307 }
0308
0309
0310 if (code == 6) {
0311 mod *= res;
0312 }
0313 return mod;
0314
0315 } else if (code == 5) {
0316 double x = paramVal;
0317 double mod = 1.0;
0318 if (x >= boundary) {
0319 mod = std::pow(high / nominal, +paramVal);
0320 } else if (x <= -boundary) {
0321 mod = std::pow(low / nominal, -paramVal);
0322 } else {
0323
0324 double x0 = boundary;
0325
0326 high /= nominal;
0327 low /= nominal;
0328
0329
0330 double logHi = std::log(high);
0331 double logLo = std::log(low);
0332 double powUp = std::exp(x0 * logHi);
0333 double powDown = std::exp(x0 * logLo);
0334 double powUpLog = high <= 0.0 ? 0.0 : powUp * logHi;
0335 double powDownLog = low <= 0.0 ? 0.0 : -powDown * logLo;
0336 double powUpLog2 = high <= 0.0 ? 0.0 : powUpLog * logHi;
0337 double powDownLog2 = low <= 0.0 ? 0.0 : -powDownLog * logLo;
0338
0339 double S0 = 0.5 * (powUp + powDown);
0340 double A0 = 0.5 * (powUp - powDown);
0341 double S1 = 0.5 * (powUpLog + powDownLog);
0342 double A1 = 0.5 * (powUpLog - powDownLog);
0343 double S2 = 0.5 * (powUpLog2 + powDownLog2);
0344 double A2 = 0.5 * (powUpLog2 - powDownLog2);
0345
0346
0347
0348 double x0Sq = x0 * x0;
0349
0350 double a = 1. / (8 * x0) * (15 * A0 - 7 * x0 * S1 + x0 * x0 * A2);
0351 double b = 1. / (8 * x0Sq) * (-24 + 24 * S0 - 9 * x0 * A1 + x0 * x0 * S2);
0352 double c = 1. / (4 * x0Sq * x0) * (-5 * A0 + 5 * x0 * S1 - x0 * x0 * A2);
0353 double d = 1. / (4 * x0Sq * x0Sq) * (12 - 12 * S0 + 7 * x0 * A1 - x0 * x0 * S2);
0354 double e = 1. / (8 * x0Sq * x0Sq * x0) * (+3 * A0 - 3 * x0 * S1 + x0 * x0 * A2);
0355 double f = 1. / (8 * x0Sq * x0Sq * x0Sq) * (-8 + 8 * S0 - 5 * x0 * A1 + x0 * x0 * S2);
0356
0357
0358 double value = 1. + x * (a + x * (b + x * (c + x * (d + x * (e + x * f)))));
0359 mod = value;
0360 }
0361 return res * (mod - 1.0);
0362 }
0363
0364 return 0.0;
0365 }
0366
0367 template <typename ParamsArray, typename DoubleArray>
0368 double flexibleInterp(unsigned int code, ParamsArray params, unsigned int n, DoubleArray low, DoubleArray high,
0369 double boundary, double nominal, int doCutoff)
0370 {
0371 double total = nominal;
0372 #if defined(__CLING__) && defined(R__HAS_CLAD)
0373 #pragma clad checkpoint loop
0374 #endif
0375 for (std::size_t i = 0; i < n; ++i) {
0376 total += flexibleInterpSingle(code, low[i], high[i], boundary, nominal, params[i], total);
0377 }
0378
0379 return doCutoff && total <= 0 ? TMath::Limits<double>::Min() : total;
0380 }
0381
0382 inline double landau(double x, double mu, double sigma)
0383 {
0384 if (sigma <= 0.)
0385 return 0.;
0386 return ROOT::Math::landau_pdf((x - mu) / sigma);
0387 }
0388
0389 inline double logNormal(double x, double k, double m0)
0390 {
0391 return ROOT::Math::lognormal_pdf(x, std::log(m0), std::abs(std::log(k)));
0392 }
0393
0394 inline double logNormalStandard(double x, double sigma, double mu)
0395 {
0396 return ROOT::Math::lognormal_pdf(x, mu, std::abs(sigma));
0397 }
0398
0399 inline double effProd(double eff, double pdf)
0400 {
0401 return eff * pdf;
0402 }
0403
0404 inline double nll(double pdf, double weight, int binnedL, int doBinOffset)
0405 {
0406 if (binnedL) {
0407
0408
0409
0410 if (std::abs(pdf) < 1e-10 && std::abs(weight) < 1e-10) {
0411 return 0.0;
0412 }
0413 if (doBinOffset) {
0414 return pdf - weight - weight * (std::log(pdf) - std::log(weight));
0415 }
0416 return pdf - weight * std::log(pdf) + TMath::LnGamma(weight + 1);
0417 } else {
0418 return -weight * std::log(pdf);
0419 }
0420 }
0421
0422 template <typename DoubleArray>
0423 double recursiveFraction(DoubleArray a, unsigned int n)
0424 {
0425 double prod = a[0];
0426
0427 for (unsigned int i = 1; i < n; ++i) {
0428 prod *= 1.0 - a[i];
0429 }
0430
0431 return prod;
0432 }
0433
0434 inline double cbShape(double m, double m0, double sigma, double alpha, double n)
0435 {
0436 double t = (m - m0) / sigma;
0437 if (alpha < 0)
0438 t = -t;
0439
0440 double absAlpha = std::abs(alpha);
0441
0442 if (t >= -absAlpha) {
0443 return std::exp(-0.5 * t * t);
0444 } else {
0445 double r = n / absAlpha;
0446 double a = std::exp(-0.5 * absAlpha * absAlpha);
0447 double b = r - absAlpha;
0448
0449 return a * std::pow(r / (b - t), n);
0450 }
0451 }
0452
0453
0454 inline double approxErf(double arg)
0455 {
0456 if (arg > 5.0)
0457 return 1.0;
0458 if (arg < -5.0)
0459 return -1.0;
0460
0461 return std::erf(arg);
0462 }
0463
0464
0465
0466
0467
0468
0469
0470
0471 inline double gaussianIntegral(double xMin, double xMax, double mean, double sigma)
0472 {
0473
0474
0475
0476 double resultScale = 0.5 * std::sqrt(TMath::TwoPi()) * sigma;
0477
0478
0479 double xscale = TMath::Sqrt2() * sigma;
0480 double scaledMin = 0.;
0481 double scaledMax = 0.;
0482 scaledMin = (xMin - mean) / xscale;
0483 scaledMax = (xMax - mean) / xscale;
0484
0485
0486
0487
0488
0489 double ecmin = std::erfc(std::abs(scaledMin));
0490 double ecmax = std::erfc(std::abs(scaledMax));
0491
0492 double cond = 0.0;
0493
0494
0495 double prd = scaledMin * scaledMax;
0496 if (prd < 0.0) {
0497 cond = 2.0 - (ecmin + ecmax);
0498 } else if (scaledMax <= 0.0) {
0499 cond = ecmax - ecmin;
0500 } else {
0501 cond = ecmin - ecmax;
0502 }
0503 return resultScale * cond;
0504 }
0505
0506 inline double bifurGaussIntegral(double xMin, double xMax, double mean, double sigmaL, double sigmaR)
0507 {
0508 const double xscaleL = TMath::Sqrt2() * sigmaL;
0509 const double xscaleR = TMath::Sqrt2() * sigmaR;
0510
0511 const double resultScale = 0.5 * std::sqrt(TMath::TwoPi());
0512
0513 if (xMax < mean) {
0514 return resultScale * (sigmaL * (std::erf((xMax - mean) / xscaleL) - std::erf((xMin - mean) / xscaleL)));
0515 } else if (xMin > mean) {
0516 return resultScale * (sigmaR * (std::erf((xMax - mean) / xscaleR) - std::erf((xMin - mean) / xscaleR)));
0517 } else {
0518 return resultScale * (sigmaR * std::erf((xMax - mean) / xscaleR) - sigmaL * std::erf((xMin - mean) / xscaleL));
0519 }
0520 }
0521
0522 inline double exponentialIntegral(double xMin, double xMax, double constant)
0523 {
0524 if (constant == 0.0) {
0525 return xMax - xMin;
0526 }
0527
0528 return (std::exp(constant * xMax) - std::exp(constant * xMin)) / constant;
0529 }
0530
0531
0532 template <bool pdfMode = false, typename DoubleArray>
0533 double polynomialIntegral(DoubleArray coeffs, int nCoeffs, int lowestOrder, double xMin, double xMax)
0534 {
0535 int denom = lowestOrder + nCoeffs;
0536 double min = coeffs[nCoeffs - 1] / double(denom);
0537 double max = coeffs[nCoeffs - 1] / double(denom);
0538
0539 for (int i = nCoeffs - 2; i >= 0; i--) {
0540 denom--;
0541 min = (coeffs[i] / double(denom)) + xMin * min;
0542 max = (coeffs[i] / double(denom)) + xMax * max;
0543 }
0544
0545 max = max * std::pow(xMax, 1 + lowestOrder);
0546 min = min * std::pow(xMin, 1 + lowestOrder);
0547
0548 return max - min + (pdfMode && lowestOrder > 0.0 ? xMax - xMin : 0.0);
0549 }
0550
0551
0552 inline double fast_fma(double x, double y, double z) noexcept
0553 {
0554 #if defined(FP_FAST_FMA)
0555 return std::fma(x, y, z);
0556 #else
0557
0558 #if defined(__clang__)
0559 #pragma STDC FP_CONTRACT ON
0560 #endif
0561 return (x * y) + z;
0562 #endif
0563 }
0564
0565 template <typename DoubleArray>
0566 double
0567 chebychevIntegral(DoubleArray coeffs, unsigned int nCoeffs, double xMin, double xMax, double xMinFull, double xMaxFull)
0568 {
0569 const double halfrange = .5 * (xMax - xMin);
0570 const double mid = .5 * (xMax + xMin);
0571
0572
0573 const double b = (xMaxFull - mid) / halfrange;
0574 const double a = (xMinFull - mid) / halfrange;
0575
0576
0577
0578
0579 double sum = b - a;
0580
0581 const unsigned int iend = nCoeffs;
0582 if (iend > 0) {
0583 {
0584
0585 const double c = coeffs[0];
0586 sum = fast_fma(0.5 * (b + a) * (b - a), c, sum);
0587 }
0588 if (1 < iend) {
0589 double bcurr = b;
0590 double btwox = 2 * b;
0591 double blast = 1;
0592
0593 double acurr = a;
0594 double atwox = 2 * a;
0595 double alast = 1;
0596
0597 double newval = atwox * acurr - alast;
0598 alast = acurr;
0599 acurr = newval;
0600
0601 newval = btwox * bcurr - blast;
0602 blast = bcurr;
0603 bcurr = newval;
0604 double nminus1 = 1.;
0605 for (unsigned int i = 1; iend != i; ++i) {
0606
0607 const double c = coeffs[i];
0608 const double term2 = (blast - alast) / nminus1;
0609
0610 newval = atwox * acurr - alast;
0611 alast = acurr;
0612 acurr = newval;
0613
0614 newval = btwox * bcurr - blast;
0615 blast = bcurr;
0616 bcurr = newval;
0617
0618 ++nminus1;
0619 const double term1 = (bcurr - acurr) / (nminus1 + 1.);
0620 const double intTn = 0.5 * (term1 - term2);
0621 sum = fast_fma(intTn, c, sum);
0622 }
0623 }
0624 }
0625
0626
0627
0628 return halfrange * sum;
0629 }
0630
0631
0632 inline double
0633 poissonIntegral(int code, double mu, double x, double integrandMin, double integrandMax, unsigned int protectNegative)
0634 {
0635 if (protectNegative && mu < 0.0) {
0636 return std::exp(-2.0 * mu);
0637 }
0638
0639 if (code == 1) {
0640
0641
0642 integrandMin = std::max(0., integrandMin);
0643
0644 if (integrandMax < 0. || integrandMax < integrandMin) {
0645 return 0;
0646 }
0647 const double delta = 100.0 * std::sqrt(mu);
0648
0649
0650
0651 if (integrandMin < std::max(mu - delta, 0.0) && integrandMax > mu + delta) {
0652 return 1.;
0653 }
0654
0655
0656 const unsigned int ixMin = integrandMin;
0657 const unsigned int ixMax = std::min(integrandMax + 1, (double)std::numeric_limits<unsigned int>::max());
0658
0659
0660 if (ixMin == 0) {
0661 return ROOT::Math::inc_gamma_c(ixMax, mu);
0662 } else {
0663
0664 if (ixMin <= mu) {
0665 return ROOT::Math::inc_gamma_c(ixMax, mu) - ROOT::Math::inc_gamma_c(ixMin, mu);
0666 } else {
0667
0668 return ROOT::Math::inc_gamma(ixMin, mu) - ROOT::Math::inc_gamma(ixMax, mu);
0669 }
0670 }
0671 }
0672
0673
0674
0675 const double ix = 1 + x;
0676
0677 return ROOT::Math::inc_gamma(ix, integrandMax) - ROOT::Math::inc_gamma(ix, integrandMin);
0678 }
0679
0680 inline double logNormalIntegral(double xMin, double xMax, double m0, double k)
0681 {
0682 const double root2 = std::sqrt(2.);
0683
0684 double ln_k = std::abs(std::log(k));
0685 double ret = 0.5 * (std::erf(std::log(xMax / m0) / (root2 * ln_k)) - std::erf(std::log(xMin / m0) / (root2 * ln_k)));
0686
0687 return ret;
0688 }
0689
0690 inline double logNormalIntegralStandard(double xMin, double xMax, double mu, double sigma)
0691 {
0692 const double root2 = std::sqrt(2.);
0693
0694 double ln_k = std::abs(sigma);
0695 double ret =
0696 0.5 * (std::erf((std::log(xMax) - mu) / (root2 * ln_k)) - std::erf((std::log(xMin) - mu) / (root2 * ln_k)));
0697
0698 return ret;
0699 }
0700
0701 inline double cbShapeIntegral(double mMin, double mMax, double m0, double sigma, double alpha, double n)
0702 {
0703 const double sqrtPiOver2 = 1.2533141373;
0704 const double sqrt2 = 1.4142135624;
0705
0706 double result = 0.0;
0707 bool useLog = false;
0708
0709 if (std::abs(n - 1.0) < 1.0e-05)
0710 useLog = true;
0711
0712 double sig = std::abs(sigma);
0713
0714 double tmin = (mMin - m0) / sig;
0715 double tmax = (mMax - m0) / sig;
0716
0717 if (alpha < 0) {
0718 double tmp = tmin;
0719 tmin = -tmax;
0720 tmax = -tmp;
0721 }
0722
0723 double absAlpha = std::abs(alpha);
0724
0725 if (tmin >= -absAlpha) {
0726 result += sig * sqrtPiOver2 * (approxErf(tmax / sqrt2) - approxErf(tmin / sqrt2));
0727 } else if (tmax <= -absAlpha) {
0728 double r = n / absAlpha;
0729 double a = r * std::exp(-0.5 * absAlpha * absAlpha);
0730 double b = r - absAlpha;
0731
0732 if (useLog) {
0733 double log_b_tmin = std::log(b - tmin);
0734 double log_b_tmax = std::log(b - tmax);
0735 result += a * std::pow(r, n - 1) * sig *
0736 (log_b_tmin - log_b_tmax + 0.5 * (1.0 - n) * (log_b_tmin * log_b_tmin - log_b_tmax * log_b_tmax));
0737 } else {
0738 result += a * sig / (1.0 - n) * (std::pow(r / (b - tmin), n - 1.0) - std::pow(r / (b - tmax), n - 1.0));
0739 }
0740 } else {
0741 double r = n / absAlpha;
0742 double a = r * std::exp(-0.5 * absAlpha * absAlpha);
0743 double b = r - absAlpha;
0744
0745 double term1 = 0.0;
0746 if (useLog) {
0747 double log_b_tmin = std::log(b - tmin);
0748 double log_r = std::log(r);
0749 term1 = a * std::pow(r, n - 1) * sig *
0750 (log_b_tmin - log_r + 0.5 * (1.0 - n) * (log_b_tmin * log_b_tmin - log_r * log_r));
0751 } else {
0752 term1 = a * sig / (1.0 - n) * (std::pow(r / (b - tmin), n - 1.0) - 1.0);
0753 }
0754
0755 double term2 = sig * sqrtPiOver2 * (approxErf(tmax / sqrt2) - approxErf(-absAlpha / sqrt2));
0756
0757 result += term1 + term2;
0758 }
0759
0760 if (result == 0)
0761 return 1.E-300;
0762 return result;
0763 }
0764
0765 template <typename DoubleArray>
0766 double bernsteinIntegral(double xlo, double xhi, double xmin, double xmax, DoubleArray coefs, int nCoefs)
0767 {
0768 double xloScaled = (xlo - xmin) / (xmax - xmin);
0769 double xhiScaled = (xhi - xmin) / (xmax - xmin);
0770
0771 int degree = nCoefs - 1;
0772 double norm = 0.;
0773
0774 for (int i = 0; i <= degree; ++i) {
0775
0776
0777
0778 double temp = 0.;
0779 for (int j = i; j <= degree; ++j) {
0780 double binCoefs = binomial(degree, j) * binomial(j, i);
0781 double oneOverJPlusOne = 1. / (j + 1.);
0782 double powDiff = std::pow(xhiScaled, j + 1.) - std::pow(xloScaled, j + 1.);
0783 temp += std::pow(-1., j - i) * binCoefs * powDiff * oneOverJPlusOne;
0784 }
0785 temp *= coefs[i];
0786 norm += temp;
0787 }
0788
0789 return norm * (xmax - xmin);
0790 }
0791
0792 template <typename XArray, typename MuArray, typename CovArray>
0793 double multiVarGaussian(int n, XArray x, MuArray mu, CovArray covI)
0794 {
0795 double result = 0.0;
0796
0797
0798 for (int i = 0; i < n; ++i) {
0799 for (int j = 0; j < n; ++j) {
0800 result += (x[i] - mu[i]) * covI[i * n + j] * (x[j] - mu[j]);
0801 }
0802 }
0803 return std::exp(-0.5 * result);
0804 }
0805
0806
0807
0808
0809 template <typename DoubleArray>
0810 double stepFunctionIntegral(double xmin, double xmax, std::size_t nBins, DoubleArray boundaries, DoubleArray coefs)
0811 {
0812 double out = 0.0;
0813 for (std::size_t i = 0; i < nBins; ++i) {
0814 double a = boundaries[i];
0815 double b = boundaries[i + 1];
0816 out += coefs[i] * std::max(0.0, std::min(b, xmax) - std::max(a, xmin));
0817 }
0818 return out;
0819 }
0820
0821 }
0822
0823 namespace clad::custom_derivatives {
0824 namespace RooFit::Detail::MathFuncs {
0825
0826
0827
0828
0829
0830
0831 template <class... Types>
0832 void binNumber_pullback(Types...)
0833 {
0834 }
0835
0836 }
0837 }
0838
0839 #endif