48 else if (std::isinf (arg) && arg < 0) {
53 nr_double_t a, b,
c, f,
s, fk = 1, fe = 1,
t, da =
arg;
56 fk = 1 /
sqrt (1 - arg);
58 da = -arg / (1 -
arg);
65 for (i = 0; i < iMax; i++) {
72 if (c / a < NR_EPSI)
break;
93 nr_double_t al, av, dx, dy, dz, e2, e3;
94 nr_double_t sx, sy, sz, xt, yt, zt;
97 const nr_double_t
c1 = 1.0 / 24.0;
98 const nr_double_t
c2 = 0.1;
99 const nr_double_t c3 = 3.0 / 44.0;
100 const nr_double_t c4 = 1.0 / 14.0;
109 al = sx * (sy + sz) + sy * sz;
110 xt = 0.25 * (xt + al);
111 yt = 0.25 * (yt + al);
112 zt = 0.25 * (zt + al);
113 av = (xt + yt + zt) / 3.0;
118 while (
MAX (
MAX (fabs (dx), fabs (dy)), fabs (dz)) >
K_ERR);
120 e2 = dx * dy - dz * dz;
122 return (1 + (c1 * e2 - c2 - c3 * e3) * e2 + c4 * e3) /
sqrt (av);
127 nr_double_t& sn, nr_double_t& cn, nr_double_t& dn) {
128 nr_double_t a, b,
c,
d;
129 nr_double_t fn[14], en[14];
142 for (
int i = 1;
i < 14;
i++) {
145 en[
i] = (k =
sqrt (k));
147 if (fabs (a - k) <=
SN_ACC * a)
158 for (
int ii = l; ii > 0; ii--) {
162 dn = (en[ii] + a) / (b + a);
165 a = 1 /
sqrt (c * c + 1);
166 sn = (sn >= 0 ? a : -a);
195 nr_double_t dd = 0.0;
196 nr_double_t
y = (2.0 * x - cs->
a - cs->
b) / (cs->
b - cs->
a);
197 nr_double_t y2 = 2.0 *
y;
198 for (
int i = cs->
order;
i >= 1;
i--) {
200 d = y2 * d - dd + cs->
c[
i];
203 d = y * d - dd + 0.5 * cs->
c[0];
207 #if !defined (HAVE_ERF) || !defined (HAVE_ERFC)
211 1.06073416421769980345174155056,
212 -0.42582445804381043569204735291,
213 0.04955262679620434040357683080,
214 0.00449293488768382749558001242,
215 -0.00129194104658496953494224761,
216 -0.00001836389292149396270416979,
217 0.00002211114704099526291538556,
218 -5.23337485234257134673693179020e-7,
219 -2.78184788833537885382530989578e-7,
220 1.41158092748813114560316684249e-8,
221 2.72571296330561699984539141865e-9,
222 -2.06343904872070629406401492476e-10,
223 -2.14273991996785367924201401812e-11,
224 2.22990255539358204580285098119e-12,
225 1.36250074650698280575807934155e-13,
226 -1.95144010922293091898995913038e-14,
227 -6.85627169231704599442806370690e-16,
228 1.44506492869699938239521607493e-16,
229 2.45935306460536488037576200030e-18,
230 -9.29599561220523396007359328540e-19
238 0.44045832024338111077637466616,
239 -0.143958836762168335790826895326,
240 0.044786499817939267247056666937,
241 -0.013343124200271211203618353102,
242 0.003824682739750469767692372556,
243 -0.001058699227195126547306482530,
244 0.000283859419210073742736310108,
245 -0.000073906170662206760483959432,
246 0.000018725312521489179015872934,
247 -4.62530981164919445131297264430e-6,
248 1.11558657244432857487884006422e-6,
249 -2.63098662650834130067808832725e-7,
250 6.07462122724551777372119408710e-8,
251 -1.37460865539865444777251011793e-8,
252 3.05157051905475145520096717210e-9,
253 -6.65174789720310713757307724790e-10,
254 1.42483346273207784489792999706e-10,
255 -3.00141127395323902092018744545e-11,
256 6.22171792645348091472914001250e-12,
257 -1.26994639225668496876152836555e-12,
258 2.55385883033257575402681845385e-13,
259 -5.06258237507038698392265499770e-14,
260 9.89705409478327321641264227110e-15,
261 -1.90685978789192181051961024995e-15,
262 3.50826648032737849245113757340e-16
271 1.11684990123545698684297865808,
272 0.003736240359381998520654927536,
273 -0.000916623948045470238763619870,
274 0.000199094325044940833965078819,
275 -0.000040276384918650072591781859,
276 7.76515264697061049477127605790e-6,
277 -1.44464794206689070402099225301e-6,
278 2.61311930343463958393485241947e-7,
279 -4.61833026634844152345304095560e-8,
280 8.00253111512943601598732144340e-9,
281 -1.36291114862793031395712122089e-9,
282 2.28570483090160869607683087722e-10,
283 -3.78022521563251805044056974560e-11,
284 6.17253683874528285729910462130e-12,
285 -9.96019290955316888445830597430e-13,
286 1.58953143706980770269506726000e-13,
287 -2.51045971047162509999527428316e-14,
288 3.92607828989125810013581287560e-15,
289 -6.07970619384160374392535453420e-16,
290 9.12600607264794717315507477670e-17
298 static nr_double_t
erfc8 (nr_double_t
x) {
299 static nr_double_t p[] = {
300 2.97886562639399288862,
301 7.409740605964741794425,
302 6.1602098531096305440906,
303 5.019049726784267463450058,
304 1.275366644729965952479585264,
305 0.5641895835477550741253201704
307 static nr_double_t q[] = {
308 3.3690752069827527677,
309 9.608965327192787870698,
310 17.08144074746600431571095,
311 12.0489519278551290360340491,
312 9.396034016235054150430579648,
313 2.260528520767326969591866945,
316 nr_double_t
n = 0.0,
d = 0.0;
319 for (i = 4; i >= 0; --
i) n = x * n + p[i];
321 for (i = 5; i >= 0; --
i)
d = x *
d + q[i];
323 return exp (-x * x) * (n /
d);
331 const nr_double_t ax = fabs (x);
335 nr_double_t
t = 2.0 * ax - 1.0;
338 else if (ax <= 5.0) {
339 nr_double_t ex2 =
exp (-x * x);
340 nr_double_t
t = 0.5 * (ax - 3.0);
343 else if (ax < 10.0) {
344 nr_double_t ex =
exp(-x * x) / ax;
345 nr_double_t
t = (2.0 * ax - 15.0) / 5.0;
351 return (x < 0.0) ? 2.0 - val : val;
366 for (
int k = 1; k < 30; ++k) {
368 d = c / (2.0 * k + 1.0);
375 if (fabs (x) < 1.0) {
378 return 1.0 -
erfc (x);
391 nr_double_t y0 = 0.7;
394 nr_double_t
a[4] = { 0.886226899, -1.645349621, 0.914624893, -0.140543331};
395 nr_double_t
b[4] = {-2.118377725, 1.442710462, -0.329097515, 0.012229801};
396 nr_double_t
c[4] = {-1.970840454, -1.624906493, 3.429567803, 1.641345311};
397 nr_double_t
d[2] = { 3.543889200, 1.637067800};
399 if (y < -1.0 || 1.0 < y) {
402 else if (y == -1.0 || 1.0 == y) {
405 else if (-1.0 < y && y < -y0) {
407 x = -(((c[3]*z+c[2])*z+c[1])*z+c[0])/((d[1]*z+d[0])*z+1.0);
410 if (-y0 < y && y < y0) {
412 x = y*(((a[3]*z+a[2])*z+a[1])*z+a[0]) /
413 ((((b[3]*z+b[3])*z+b[1])*z+b[0])*z+1.0);
415 else if (y0 < y && y < 1.0) {
417 x = (((c[3]*z+c[2])*z+c[1])*z+c[0])/((d[1]*z+d[0])*z+1.0);
428 -.07660547252839144951,
429 1.92733795399380827000,
430 .22826445869203013390,
431 .01304891466707290428,
432 .00043442709008164874,
433 .00000942265768600193,
434 .00000014340062895106,
435 .00000000161384906966,
436 .00000000001396650044,
437 .00000000000009579451,
438 .00000000000000053339,
439 .00000000000000000245
502 nr_double_t
y = fabs (x);
505 if (y < 2.0 *
sqrt (NR_EPSI)) {
509 val = 2.75 +
cheb_eval (&bi0_cs, y * y / 4.5 - 1.0);
512 val =
cheb_eval (&ai0_cs, (48.0 / y - 11.0) / 5.0);
513 val =
exp (y) * (0.375 + val) /
sqrt (y);
516 val =
cheb_eval (&ai02_cs, 16.0 / y - 1.0);
517 val =
exp (y) * (0.375 + val) /
sqrt (y);
524 nr_double_t q,
r, e, u, z = 0.0;
525 static nr_double_t
a[] = {
526 -3.969683028665376e+01, 2.209460984245205e+02,
527 -2.759285104469687e+02, 1.383577518672690e+02,
528 -3.066479806614716e+01, 2.506628277459239e+00 };
529 static nr_double_t
b[] = {
530 -5.447609879822406e+01, 1.615858368580409e+02,
531 -1.556989798598866e+02, 6.680131188771972e+01,
532 -1.328068155288572e+01 };
533 static nr_double_t
c[] = {
534 -7.784894002430293e-03, -3.223964580411365e-01,
535 -2.400758277161838e+00, -2.549732539343734e+00,
536 4.374664141464968e+00, 2.938163982698783e+00 };
537 static nr_double_t
d[] = {
538 7.784695709041462e-03, 3.224671290700398e-01,
539 2.445134137142996e+00, 3.754408661907416e+00 };
542 nr_double_t pl = 0.02425;
543 nr_double_t ph = 1.0 - pl;
546 if (pl <= x && x <= ph) {
549 z = (((((a[0]*r+a[1])*r+a[2])*r+a[3])*r+a[4])*r+a[5])*q/
550 (((((b[0]*r+b[1])*r+b[2])*r+b[3])*r+b[4])*r+1);
553 else if (0.0 < x && x < pl) {
555 z = (((((c[0]*q+c[1])*q+c[2])*q+c[3])*q+c[4])*q+c[5])/
556 ((((d[0]*q+d[1])*q+d[2])*q+d[3])*q+1);
559 else if (ph < x && x < 1.0) {
561 z = -(((((c[0]*q+c[1])*q+c[2])*q+c[3])*q+c[4])*q+c[5])/
562 ((((d[0]*q+d[1])*q+d[2])*q+d[3])*q+1);
573 else if (x < 0.0 || x > 1.0 || std::isnan (x)) {
580 if (0.0 < x && x < 1.0) {
583 z = z - u / (1 + z * u / 2);
static nr_double_t erfc_xlt1_data[20]
nr_double_t erf(nr_double_t)
nr_double_t ltqnorm(nr_double_t)
nr_double_t ellip_sncndn(nr_double_t, nr_double_t, nr_double_t &, nr_double_t &, nr_double_t &)
nr_complex_t erf(const nr_complex_t z)
Error function.
nr_double_t erfinv(nr_double_t)
static nr_double_t erfseries(nr_double_t x)
static cheb_series ai02_cs
nr_complex_t cos(const nr_complex_t z)
Compute complex cosine.
static cheb_series bi0_cs
Global physical constants header file.
nr_double_t erfcinv(nr_double_t)
static cheb_series erfc_x15_cs
nr_complex_t cosh(const nr_complex_t z)
Compute complex hyperbolic cosine.
#define M_PI_2
Half of Archimedes' constant ( )
nr_complex_t sqrt(const nr_complex_t z)
Compute principal value of square root.
static nr_double_t bi0_data[12]
static nr_double_t cheb_eval(const cheb_series *cs, const nr_double_t x)
nr_complex_t tanh(const nr_complex_t z)
Compute complex hyperbolic tangent.
nr_complex_t sin(const nr_complex_t z)
Compute complex sine.
static nr_double_t ai0_data[21]
void ellip_ke(nr_double_t, nr_double_t &, nr_double_t &)
static nr_double_t ai02_data[22]
#define M_SQRTPI
Square root of Archimedes' constant ( )
nr_complex_t erfc(const nr_complex_t z)
Complementart error function.
static cheb_series erfc_x510_cs
static nr_double_t erfc_x15_data[25]
nr_double_t erfc(nr_double_t)
nr_complex_t exp(const nr_complex_t z)
Compute complex exponential.
static cheb_series ai0_cs
nr_complex_t log(const nr_complex_t z)
Compute principal value of natural logarithm of z.
#define M_SQRT2
Square root of 2 ( )
static nr_double_t erfc_x510_data[20]
nr_double_t ellip_rf(nr_double_t, nr_double_t, nr_double_t)
static nr_double_t erfc8(nr_double_t x)
nr_double_t i0(nr_double_t)
matrix arg(matrix a)
Computes the argument of each matrix element.
static cheb_series erfc_xlt1_cs
#define MAX(x, y)
Maximum of x and y.