54 Set (kappa, beta2, epsilonPM, epsilon);
78 static const double eu = 0.577215664901532860606;
79 static const double pi2 = 6.28318530717958647693,
80 rpi = 0.318309886183790671538,
81 pih = 1.57079632679489661923;
83 double deltaEpsilon = 0.001;
84 static const double logdeltaEpsilon = -
std::log(deltaEpsilon);
86 static const double eps = 1
e-5;
89 9.29, 2.47, 0.89, 0.36, 0.15, 0.07, 0.03, 0.02};
91 0.012, 0.03, 0.08, 0.26, 0.87, 3.83};
94 std::cerr <<
"VavilovAccurate::Set: kappa = " << kappa <<
" - out of range" << std::endl;
95 if (kappa < 0.001) kappa = 0.001;
97 if (beta2 < 0 || beta2 > 1) {
98 std::cerr <<
"VavilovAccurate::Set: beta2 = " << beta2 <<
" - out of range" << std::endl;
99 if (beta2 < 0) beta2 = -beta2;
100 if (beta2 > 1) beta2 = 1;
104 fH[5] = 1-beta2*(1-
eu)-logEpsilonPM/kappa;
107 double h4 = logEpsilonPM/kappa-(1+beta2*
eu);
109 double kappaInv = 1/kappa;
115 while (lp < 9 && kappa < xp[lp]) ++lp;
117 while (lq < 7 && kappa >= xq[lq]) ++
lq;
124 }
while (ifail == 2);
134 fH[1] = kappa*(2+beta2*
eu)+h1;
135 if(kappa >= 0.07) fH[1] += logdeltaEpsilon;
157 double d = rpi*
std::exp(kappa*(1+beta2*(eu-logKappa)));
162 for (
int k = 1; k <
n; ++k) {
165 double x1 = kappaInv*
x;
170 double xf1 = kappa*(beta2*c1-c4)-x*c2;
171 double xf2 = x*(c1 +
fT0) + kappa*(c3+beta2*c2);
172 double d1 = q*d*fOmega*
std::exp(xf1);
191 if (
fKappa < 0.02)
return;
197 if (estmedian>1.3) estmedian = 1.3;
200 for (
int i = 1; i <
fNQuant/2; ++i) {
206 double x = estmedian + (i-fNQuant/2)*(
fT1-estmedian)/(fNQuant/2-1);
220 static const double pi = 3.14159265358979323846;
226 }
else if (x <=
fT1) {
229 double cof = 2*
cos(u);
233 for (
int k = 2; k <= n+1; ++k) {
240 for (
int k = 2; k <=
n; ++k) {
245 f = 0.5*(a0-a2)+b0*
sin(u);
259 static const double pi = 3.14159265358979323846;
265 }
else if (x <=
fT1) {
268 double cof = 2*
cos(u);
272 for (
int k = 2; k <= n+1; ++k) {
279 for (
int k = 2; k <=
n; ++k) {
284 f = 0.5*(a0-a2)+b0*
sin(u);
299 static const double pi = 3.14159265358979323846;
305 }
else if (x <=
fT1) {
308 double cof = 2*
cos(u);
312 for (
int k = 2; k <= n+1; ++k) {
319 for (
int k = 2; k <=
n; ++k) {
324 f = -0.5*(a0-a2)+b0*
sin(u);
338 if (z < 0 || z > 1)
return std::numeric_limits<double>::signaling_NaN();
351 while (z >
fQuant[i]) ++i;
367 assert (x >
fT0 && x <
fT1);
377 if (x <
fT0) x = 0.5*(
fT0+x-dx);
378 else if (x >
fT1) x = 0.5*(
fT1+x-dx);
379 assert (x >
fT0 && x <
fT1);
390 if (z < 0 || z > 1)
return std::numeric_limits<double>::signaling_NaN();
405 while (z1 >
fQuant[i]) ++i;
429 assert (x >
fT0 && x <
fT1);
435 double y1 = -
Pdf (x);
439 if (x <
fT0) x = 0.5*(
fT0+x-dx);
440 else if (x >
fT1) x = 0.5*(
fT1+x-dx);
441 assert (x >
fT0 && x <
fT1);
464 return vavilov->
Pdf (x);
469 return vavilov->
Cdf_c (x);
474 return vavilov->
Cdf (x);
508 double xa, xb, fa, fb,
r;
529 bool recalcF12 =
true;
530 bool recalcFab =
true;
534 double x1=0,
x2=0,
f1=0,
f2=0, fx=0, ee=0;
539 ee = eps*(std::abs(x0)+1);
568 if(u2 == 0 || u4 == 0)
continue;
574 double cb = (x1+
x2)*u2-(x2+x0)*u1;
575 double cc = (x1-x0)*
f1-x1*(ca*x1+cb);
577 if(cb == 0)
continue;
583 x0 = -u3 + (x0+u3 >= 0 ? +1 : -1)*
std::sqrt(u4);
585 if(x0 < xa || x0 > xb)
continue;
589 r = std::abs(x0-x3) < std::abs(x0-x2) ? std::abs(x0-x3) : std::abs(x0-x2);
590 ee = eps*(std::abs(x0)+1);
601 fx = (this->*
f) (x0);
617 if (fx*ff <= 0)
break;
633 r = -0.5*std::abs(xb-xa);
635 std::cerr <<
"VavilovAccurate::Rzero: fail=" << fail <<
", f(" << a <<
")=" << (this->*
f) (a)
636 <<
", f(" << b <<
")=" << (this->*
f) (b) << std::endl;
646 static const double eu = 0.577215664901532860606;
649 return (x-0.25*x)*x-
eu;
678 if (x>-0.223172) x = -0.223172;
683 double p0 =
Pdf (x - eps);
685 double p2 =
Pdf (x + eps);
686 double y1 = 0.5*(p2-p0)/eps;
687 double y2 = (p2-2*p1+p0)/(eps*eps);
690 if (
fabs(dx) < eps) eps = 0.1*
fabs(dx);
691 }
while (
fabs(dx) > 1
E-5);
double cosint(double x)
Calculates the real part of the cosine integral Re(Ci).
virtual double GetKappa() const
Return the current value of .
Namespace for new ROOT classes and functions.
double sinint(double x)
Calculates the sine integral.
int Rzero(double a, double b, double &x0, double eps, int mxf, double(VavilovAccurate::*f)(double) const) const
double fA_cdf[MAXTERMS+1]
double G116f2(double x) const
double Cdf(double x) const
Evaluate the Vavilov cumulative probability density function.
double GetNTerms() const
Return the number of terms used in the series expansion.
double vavilov_accurate_cdf_c(double x, double kappa, double beta2)
The Vavilov complementary cumulative probability density function.
double GetEpsilonPM() const
Return the current value of .
static const double x2[5]
virtual void SetKappaBeta2(double kappa, double beta2)
Change and and recalculate coefficients if necessary.
static double p2(double t, double a, double b, double c)
double fQuant[kNquantMax]
double fB_cdf[MAXTERMS+1]
virtual double GetLambdaMax() const
Return the maximum value of for which is nonzero in the current approximation.
double GetEpsilon() const
Return the current value of .
double expint(double x)
Calculates the exponential integral.
void Set(double kappa, double beta2, double epsilonPM=5E-4, double epsilon=1E-5)
(Re)Initialize the object
static VavilovAccurate * GetInstance()
Returns a static instance of class VavilovFast.
VecExpr< UnaryOp< Fabs< T >, VecExpr< A, T, D >, T >, T, D > fabs(const VecExpr< A, T, D > &rhs)
double Quantile(double z) const
Evaluate the inverse of the Vavilov cumulative probability density function.
double Cdf_c(double x) const
Evaluate the Vavilov complementary cumulative probability density function.
virtual double Mode() const
Return the value of where the pdf is maximal.
static VavilovAccurate * fgInstance
double fA_pdf[MAXTERMS+1]
static double p1(double t, double a, double b)
VavilovAccurate(double kappa=1, double beta2=1, double epsilonPM=5E-4, double epsilon=1E-5)
Initialize an object to calculate the Vavilov distribution.
void InitQuantile() const
static const double x1[5]
double vavilov_accurate_pdf(double x, double kappa, double beta2)
The Vavilov probability density function.
Class describing a Vavilov distribution.
virtual double GetLambdaMin() const
Return the minimum value of for which is nonzero in the current approximation.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t e
double G116f1(double x) const
Namespace for new Math classes and functions.
you should not use this method at all Int_t Int_t z
virtual double GetBeta2() const
Return the current value of .
static double E1plLog(double x)
double landau_quantile(double z, double xi=1)
Inverse ( ) of the cumulative distribution function of the lower tail of the Landau distribution (lan...
double Pdf(double x) const
Evaluate the Vavilov probability density function.
double f2(const double *x)
double fLambda[kNquantMax]
double Quantile_c(double z) const
Evaluate the inverse of the complementary Vavilov cumulative probability density function.
virtual ~VavilovAccurate()
Destructor.
you should not use this method at all Int_t Int_t Double_t Double_t Double_t Int_t Double_t Double_t Double_t Double_t b
double fB_pdf[MAXTERMS+1]
double vavilov_accurate_quantile_c(double z, double kappa, double beta2)
The inverse of the complementary Vavilov cumulative probability density function. ...
double vavilov_accurate_quantile(double z, double kappa, double beta2)
The inverse of the Vavilov cumulative probability density function.
double vavilov_accurate_cdf(double x, double kappa, double beta2)
The Vavilov cumulative probability density function.
static const double x3[11]