Rosetta
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Properties Friends Macros Pages
functions.cc
Go to the documentation of this file.
1 // -*- mode:c++;tab-width:2;indent-tabs-mode:t;show-trailing-whitespace:t;rm-trailing-spaces:t -*-
2 // vi: set ts=2 noet:
3 //
4 // (c) Copyright Rosetta Commons Member Institutions.
5 // (c) This file is part of the Rosetta software suite and is made available under license.
6 // (c) The Rosetta software is developed by the contributing members of the Rosetta Commons.
7 // (c) For more information, see http://www.rosettacommons.org. Questions about this can be
8 // (c) addressed to University of Washington UW TechTransfer, email: license@u.washington.edu.
9 
10 /// @file numeric/statistics.functions.cc
11 /// @brief
12 /// @author James Thompson
13 
14 // Platform headers
16 #include <numeric/types.hh>
17 #include <utility/vector1.hh>
18 #include <utility/numbers.hh>
19 
20 #include <cmath>
21 #include <cfloat>
22 #include <limits>
23 
24 
25 // defines for error functions
26 typedef std::complex<double> cmplx;
27 #define C(a,b) cmplx(a,b)
28 #define cexp(z) std::exp(z)
29 #define creal(z) std::real(z)
30 #define cimag(z) std::imag(z)
31 #define cpolar(r,t) std::polar(r,t)
32 #define Inf std::numeric_limits<double>::infinity()
33 #define NaN std::numeric_limits<double>::quiet_NaN()
34 
35 namespace numeric {
36 namespace statistics {
37 
40  utility::vector1< numeric::Real > const & posterior
41 ) {
42  assert( prior.size() == posterior.size() );
43 
44  // p = prior distribution, q = posterior distribution
45  // relative entropy = sum( p(i) * log( p(i) / q(i) )
46  Real div( 0.0 );
48  for ( iter p = prior.begin(), p_end = prior.end(),
49  q = posterior.begin(), q_end = posterior.end();
50  p != p_end && q != q_end; ++p, ++q
51  ) {
52  div += *p * std::log( *p / *q );
53  }
54 
55  return div;
56 }
57 
61 
62 {
63  assert( vec1.size() == vec2.size() );
64  numeric::Real m1 = mean( vec1.begin(),vec1.end(),0.0 );
65  numeric::Real m2 = mean( vec2.begin(),vec2.end(),0.0 );
66  numeric::Real sd1 = std_dev_with_provided_mean( vec1.begin(), vec1.end(), m1 );
67  numeric::Real sd2 = std_dev_with_provided_mean( vec2.begin(), vec2.end(), m2 );
68  return corrcoef_with_provided_mean_and_std_dev( vec1, m1, sd1, vec2, m2, sd2 );
69 }
70 
73  numeric::Real m1,
74  numeric::Real sd1,
76  numeric::Real m2,
77  numeric::Real sd2
78 ) {
79  numeric::Real cov(0);
80  Size n(0);
82  for ( iter v1 = vec1.begin(), v2 = vec2.begin(), v1_end = vec1.end(), v2_end = vec2.end();
83  v1 != v1_end && v2 != v2_end;
84  ++v1, ++v2, n++ ) {
85  cov += ( *v1 - m1 ) * ( *v2 - m2 );
86  }
87  cov /= ( n - 1 );
88  return cov / ( sd1 * sd2 );
89 }
90 
94 ) {
95  assert( vec1.size() == vec2.size() );
96  numeric::Real m1 = mean( vec1.begin(),vec1.end(),0.0 );
97  numeric::Real m2 = mean( vec2.begin(),vec2.end(),0.0 );
98  return cov_with_provided_mean( vec1, m1, vec2, m2 );
99 }
100 
103  numeric::Real m1,
105  numeric::Real m2
106 ) {
107  numeric::Real cov(0);
108  Size n(0);
110  for ( iter v1 = vec1.begin(), v2 = vec2.begin(), v1_end = vec1.end(), v2_end = vec2.end();
111  v1 != v1_end && v2 != v2_end;
112  ++v1, ++v2, n++ ) {
113  cov += ( *v1 - m1 ) * ( *v2 - m2 );
114  }
115  cov /= ( n - 1 );
116  return cov;
117 }
118 
119 
120 //fpd error function with complex argmuments
121 // code from http://ab-initio.mit.edu/Faddeeva
122 // Auxiliary routines to compute other special functions based on w(z)
123 
124 
125 // compute erfcx(z) = exp(z^2) erfz(z)
126 cmplx errfcx(cmplx z, double relerr) {
127  return w(C(-cimag(z), creal(z)), relerr);
128 }
129 
130 // compute the error function erf(x)
131 double errf(double x)
132 {
133  double mx2 = -x*x;
134  if ( mx2 < -750 ) { // underflow
135  return (x >= 0 ? 1.0 : -1.0);
136  }
137 
138  if ( x >= 0 ) {
139  if ( x < 8e-2 ) goto taylor;
140  return 1.0 - exp(mx2) * errfcx(x);
141  } else { // x < 0
142  if ( x > -8e-2 ) goto taylor;
143  return exp(mx2) * errfcx(-x) - 1.0;
144  }
145 
146  // Use Taylor series for small |x|, to avoid cancellation inaccuracy
147  // erf(x) = 2/sqrt(pi) * x * (1 - x^2/3 + x^4/10 - x^6/42 + x^8/216 + ...)
148  taylor:
149  return x * (1.1283791670955125739
150  + mx2 * (0.37612638903183752464
151  + mx2 * (0.11283791670955125739
152  + mx2 * (0.026866170645131251760
153  + mx2 * 0.0052239776254421878422))));
154 }
155 
156 // compute the error function erf(z)
157 cmplx errf(cmplx z, double relerr)
158 {
159  double x = creal(z), y = cimag(z);
160 
161  if ( y == 0 ) {
162  return C(errf(x), y); // preserve sign of 0
163  }
164  if ( x == 0 ) { // handle separately for speed & handling of y = Inf or NaN
165  //handle y -> Inf limit manually, since
166  // exp(y^2) -> Inf but Im[w(y)] -> 0, so
167  // IEEE will give us a NaN when it should be Inf
168  return C(x, y*y > 720 ? (y > 0 ? Inf : -Inf): exp(y*y) * w_im(y));
169  }
170 
171  double mRe_z2 = (y - x) * (x + y); // Re(-z^2), being careful of overflow
172  double mIm_z2 = -2*x*y; // Im(-z^2)
173  if ( mRe_z2 < -750 ) { // underflow
174  return (x >= 0 ? 1.0 : -1.0);
175  }
176 
177  // Handle positive and negative x via different formulas,
178  // using the mirror symmetries of w, to avoid overflow/underflow
179  // problems from multiplying exponentially large and small quantities.
180  if ( x >= 0 ) {
181  if ( x < 8e-2 ) {
182  if ( fabs(y) < 1e-2 ) {
183  goto taylor;
184  } else if ( fabs(mIm_z2) < 5e-3 && x < 5e-3 ) {
185  goto taylor_erfi;
186  }
187  }
188  // don't use complex exp function, since that will produce spurious NaN
189  // values when multiplying w in an overflow situation.
190  return 1.0 - exp(mRe_z2) *
191  (C(cos(mIm_z2), sin(mIm_z2)) * w(C(-y,x), relerr));
192  } else { // x < 0
193  if ( x > -8e-2 ) { // duplicate from above to avoid fabs(x) call
194  if ( fabs(y) < 1e-2 ) {
195  goto taylor;
196  } else if ( fabs(mIm_z2) < 5e-3 && x > -5e-3 ) {
197  goto taylor_erfi;
198  }
199  } else if ( utility::is_nan(x) ) {
200  return C(NaN, y == 0 ? 0 : NaN);
201  }
202  // don't use complex exp function, since that will produce spurious NaN
203  // values when multiplying w in an overflow situation. */
204  return exp(mRe_z2) *
205  (C(cos(mIm_z2), sin(mIm_z2))
206  * w(C(y,-x), relerr)) - 1.0;
207  }
208 
209  // Use Taylor series for small |z|, to avoid cancellation inaccuracy
210  // erf(z) = 2/sqrt(pi) * z * (1 - z^2/3 + z^4/10 - z^6/42 + z^8/216 + ...)
211  taylor:
212  {
213  cmplx mz2 = C(mRe_z2, mIm_z2); // -z^2
214  return z * (1.1283791670955125739
215  + mz2 * (0.37612638903183752464
216  + mz2 * (0.11283791670955125739
217  + mz2 * (0.026866170645131251760
218  + mz2 * 0.0052239776254421878422))));
219  }
220 
221  taylor_erfi:
222  {
223  double x2 = x*x, y2 = y*y;
224  double expy2 = exp(y2);
225  return C
226  (expy2 * x * (1.1283791670955125739
227  - x2 * (0.37612638903183752464
228  + 0.75225277806367504925*y2)
229  + x2*x2 * (0.11283791670955125739
230  + y2 * (0.45135166683820502956
231  + 0.15045055561273500986*y2))),
232  expy2 * (w_im(y)
233  - x2*y * (1.1283791670955125739
234  - x2 * (0.56418958354775628695
235  + 0.37612638903183752464*y2))));
236  }
237 }
238 
239 // erfi(z) = -i erf(iz)
240 cmplx errfi(cmplx z, double relerr) {
241  cmplx e = errf(C(-cimag(z),creal(z)), relerr);
242  return C(cimag(e), -creal(e));
243 }
244 
245 // erfi(x) = -i erf(ix)
246 double errfi(double x) {
247  return x*x > 720 ? (x > 0 ? Inf : -Inf)
248  : exp(x*x) * w_im(x);
249 }
250 
251 // erfc(x) = 1 - erf(x)
252 double errfc(double x) {
253  if ( x*x > 750 ) { // underflow
254  return (x >= 0 ? 0.0 : 2.0);
255  }
256  return x >= 0 ? exp(-x*x) * errfcx(x) : 2. - exp(-x*x) * errfcx(-x);
257 }
258 
259 // erfc(z) = 1 - erf(z)
260 cmplx errfc(cmplx z, double relerr) {
261  double x = creal(z), y = cimag(z);
262 
263  if ( x == 0. ) {
264  return C(1,
265  // handle y -> Inf limit manually, since
266  // exp(y^2) -> Inf but Im[w(y)] -> 0, so
267  // IEEE will give us a NaN when it should be Inf
268  y*y > 720 ? (y > 0 ? -Inf : Inf)
269  : -exp(y*y) * w_im(y));
270  }
271  if ( y == 0. ) {
272  if ( x*x > 750 ) { // underflow
273  return C(x >= 0 ? 0.0 : 2.0,
274  -y); // preserve sign of 0
275  }
276  return C(x >= 0 ? exp(-x*x) * errfcx(x)
277  : 2. - exp(-x*x) * errfcx(-x),
278  -y); // preserve sign of zero
279  }
280 
281  double mRe_z2 = (y - x) * (x + y); // Re(-z^2), being careful of overflow
282  double mIm_z2 = -2*x*y; // Im(-z^2)
283  if ( mRe_z2 < -750 ) { // underflow
284  return (x >= 0 ? 0.0 : 2.0);
285  }
286 
287  if ( x >= 0 ) {
288  return cexp(C(mRe_z2, mIm_z2))
289  * w(C(-y,x), relerr);
290  } else {
291  return 2.0 - cexp(C(mRe_z2, mIm_z2))
292  * w(C(y,-x), relerr);
293  }
294 }
295 
296 // compute Dawson(x) = sqrt(pi)/2 * exp(-x^2) * erfi(x)
297 double Dawson(double x) {
298  const double spi2 = 0.8862269254527580136490837416705725913990; // sqrt(pi)/2
299  return spi2 * w_im(x);
300 }
301 
302 // compute Dawson(z) = sqrt(pi)/2 * exp(-z^2) * erfi(z)
303 cmplx Dawson(cmplx z, double relerr) {
304  const double spi2 = 0.8862269254527580136490837416705725913990; // sqrt(pi)/2
305  double x = creal(z), y = cimag(z);
306 
307  // handle axes separately for speed & proper handling of x or y = Inf or NaN
308  if ( y == 0 ) {
309  return C(spi2 * w_im(x),
310  -y); // preserve sign of 0
311  }
312  if ( x == 0 ) {
313  double y2 = y*y;
314  if ( y2 < 2.5e-5 ) { // Taylor expansion
315  return C(x, // preserve sign of 0
316  y * (1.
317  + y2 * (0.6666666666666666666666666666666666666667
318  + y2 * 0.26666666666666666666666666666666666667)));
319  }
320  return C(x, // preserve sign of 0
321  spi2 * (y >= 0
322  ? exp(y2) - errfcx(y) : errfcx(-y) - exp(y2)));
323  }
324 
325  double mRe_z2 = (y - x) * (x + y); // Re(-z^2), being careful of overflow
326  double mIm_z2 = -2*x*y; // Im(-z^2)
327  cmplx mz2 = C(mRe_z2, mIm_z2); // -z^2
328 
329  // Handle positive and negative x via different formulas,
330  // using the mirror symmetries of w, to avoid overflow/underflow
331  // problems from multiplying exponentially large and small quantities.
332  if ( y >= 0 ) {
333  if ( y < 5e-3 ) {
334  if ( fabs(x) < 5e-3 ) {
335  goto taylor;
336  } else if ( fabs(mIm_z2) < 5e-3 ) {
337  goto taylor_realaxis;
338  }
339  }
340  cmplx res = cexp(mz2) - w(z, relerr);
341  return spi2 * C(-cimag(res), creal(res));
342  } else { // y < 0
343  if ( y > -5e-3 ) { // duplicate from above to avoid fabs(x) call
344  if ( fabs(x) < 5e-3 ) {
345  goto taylor;
346  } else if ( fabs(mIm_z2) < 5e-3 ) {
347  goto taylor_realaxis;
348  }
349  } else if ( utility::is_nan(y) ) {
350  return C(x == 0 ? 0 : NaN, NaN);
351  }
352  cmplx res = w(-z, relerr) - cexp(mz2);
353  return spi2 * C(-cimag(res), creal(res));
354  }
355 
356  // Use Taylor series for small |z|, to avoid cancellation inaccuracy
357  // dawson(z) = z - 2/3 z^3 + 4/15 z^5 + ...
358  taylor:
359  return z * (1.
360  + mz2 * (0.6666666666666666666666666666666666666667
361  + mz2 * 0.2666666666666666666666666666666666666667));
362 
363  // for small |y| and small |xy|,
364  // use Taylor series to avoid cancellation inaccuracy:
365  // dawson(x + iy)
366  // = D + y^2 (D + x - 2Dx^2)
367  // + y^4 (D/2 + 5x/6 - 2Dx^2 - x^3/3 + 2Dx^4/3)
368  // + iy [ (1-2Dx) + 2/3 y^2 (1 - 3Dx - x^2 + 2Dx^3)
369  // + y^4/15 (4 - 15Dx - 9x^2 + 20Dx^3 + 2x^4 - 4Dx^5) ] + ...
370  // where D = dawson(x)
371  //
372  // However, for large |x|, 2Dx -> 1 which gives cancellation problems in
373  // this series (many of the leading terms cancel). So, for large |x|,
374  // we need to substitute a continued-fraction expansion for D.
375  //
376  // dawson(x) = 0.5 / (x-0.5/(x-1/(x-1.5/(x-2/(x-2.5/(x...))))))
377  //
378  // The 6 terms shown here seems to be the minimum needed to be
379  // accurate as soon as the simpler Taylor expansion above starts
380  // breaking down. Using this 6-term expansion, factoring out the
381  // denominator, and simplifying with Maple, we obtain:
382  //
383  // Re dawson(x + iy) * (-15 + 90x^2 - 60x^4 + 8x^6) / x
384  // = 33 - 28x^2 + 4x^4 + y^2 (18 - 4x^2) + 4 y^4
385  // Im dawson(x + iy) * (-15 + 90x^2 - 60x^4 + 8x^6) / y
386  // = -15 + 24x^2 - 4x^4 + 2/3 y^2 (6x^2 - 15) - 4 y^4
387  //
388  // Finally, for |x| > 5e7, we can use a simpler 1-term continued-fraction
389  // expansion for the real part, and a 2-term expansion for the imaginary
390  // part. (This avoids overflow problems for huge |x|.) This yields:
391  //
392  // Re dawson(x + iy) = [1 + y^2 (1 + y^2/2 - (xy)^2/3)] / (2x)
393  // Im dawson(x + iy) = y [ -1 - 2/3 y^2 + y^4/15 (2x^2 - 4) ] / (2x^2 - 1)
394  taylor_realaxis:
395  {
396  double x2 = x*x;
397  if ( x2 > 1600 ) { // |x| > 40
398  double y2 = y*y;
399  if ( x2 > 25e14 ) { // |x| > 5e7
400  double xy2 = (x*y)*(x*y);
401  return C((0.5 + y2 * (0.5 + 0.25*y2
402  - 0.16666666666666666667*xy2)) / x,
403  y * (-1 + y2 * (-0.66666666666666666667
404  + 0.13333333333333333333*xy2
405  - 0.26666666666666666667*y2))
406  / (2*x2 - 1));
407  }
408  return (1. / (-15 + x2*(90 + x2*(-60 + 8*x2)))) *
409  C(x * (33 + x2 * (-28 + 4*x2)
410  + y2 * (18 - 4*x2 + 4*y2)),
411  y * (-15 + x2 * (24 - 4*x2)
412  + y2 * (4*x2 - 10 - 4*y2)));
413  } else {
414  double D = spi2 * w_im(x);
415  double y2 = y*y;
416  return C
417  (D + y2 * (D + x - 2*D*x2)
418  + y2*y2 * (D * (0.5 - x2 * (2 - 0.66666666666666666667*x2))
419  + x * (0.83333333333333333333
420  - 0.33333333333333333333 * x2)),
421  y * (1 - 2*D*x
422  + y2 * 0.66666666666666666667 * (1 - x2 - D*x * (3 - 2*x2))
423  + y2*y2 * (0.26666666666666666667 -
424  x2 * (0.6 - 0.13333333333333333333 * x2)
425  - D*x * (1 - x2 * (1.3333333333333333333
426  - 0.26666666666666666667 * x2)))));
427  }
428  }
429 }
430 
431 /////////////////////////////////////////////////////////////////////////
432 
433 // return sinc(x) = sin(x)/x, given both x and sin(x)
434 // [since we only use this in cases where sin(x) has already been computed]
435 static inline double sinc(double x, double sinx) {
436  return fabs(x) < 1e-4 ? 1 - (0.1666666666666666666667)*x*x : sinx / x;
437 }
438 
439 // sinh(x) via Taylor series, accurate to machine precision for |x| < 1e-2
440 static inline double sinh_taylor(double x) {
441  return x * (1 + (x*x) * (0.1666666666666666666667
442  + 0.00833333333333333333333 * (x*x)));
443 }
444 
445 static inline double sqr(double x) { return x*x; }
446 
447 // precomputed table of expa2n2[n-1] = exp(-a2*n*n)
448 // for double-precision a2 = 0.26865... in w, below.
449 static const double expa2n2[] = {
450 7.64405281671221563e-01, 3.41424527166548425e-01, 8.91072646929412548e-02, 1.35887299055460086e-02,
451 1.21085455253437481e-03, 6.30452613933449404e-05, 1.91805156577114683e-06, 3.40969447714832381e-08,
452 3.54175089099469393e-10, 2.14965079583260682e-12, 7.62368911833724354e-15, 1.57982797110681093e-17,
453 1.91294189103582677e-20, 1.35344656764205340e-23, 5.59535712428588720e-27, 1.35164257972401769e-30,
454 1.90784582843501167e-34, 1.57351920291442930e-38, 7.58312432328032845e-43, 2.13536275438697082e-47,
455 3.51352063787195769e-52, 3.37800830266396920e-57, 1.89769439468301000e-62, 6.22929926072668851e-68,
456 1.19481172006938722e-73, 1.33908181133005953e-79, 8.76924303483223939e-86, 3.35555576166254986e-92,
457 7.50264110688173024e-99, 9.80192200745410268e-106, 7.48265412822268959e-113, 3.33770122566809425e-120,
458 8.69934598159861140e-128, 1.32486951484088852e-135, 1.17898144201315253e-143, 6.13039120236180012e-152,
459 1.86258785950822098e-160, 3.30668408201432783e-169, 3.43017280887946235e-178, 2.07915397775808219e-187,
460 7.36384545323984966e-197, 1.52394760394085741e-206, 1.84281935046532100e-216, 1.30209553802992923e-226,
461 5.37588903521080531e-237, 1.29689584599763145e-247, 1.82813078022866562e-258, 1.50576355348684241e-269,
462 7.24692320799294194e-281, 2.03797051314726829e-292, 3.34880215927873807e-304,
463 0.0 // underflow (also prevents reads past array end, below)
464 };
465 
466 static const unsigned int expa2n2_length = 52;
467 
468 /////////////////////////////////////////////////////////////////////////
469 
470 cmplx w(cmplx z, double relerr) {
471  if ( creal(z) == 0.0 ) {
472  return C(errfcx(cimag(z)), creal(z)); // give correct sign of 0 in cimag(w)
473  } else if ( cimag(z) == 0 ) {
474  return C(exp(-sqr(creal(z))), w_im(creal(z)));
475  }
476 
477  double a, a2, c;
478  if ( relerr <= DBL_EPSILON ) {
479  relerr = DBL_EPSILON;
480  a = 0.518321480430085929872; // pi / sqrt(-log(eps*0.5))
481  c = 0.329973702884629072537; // (2/pi) * a;
482  a2 = 0.268657157075235951582; // a^2
483  } else {
484  const double pi = 3.14159265358979323846264338327950288419716939937510582;
485  if ( relerr > 0.1 ) relerr = 0.1; // not sensible to compute < 1 digit
486  a = pi / sqrt(-log(relerr*0.5));
487  c = (2/pi)*a;
488  a2 = a*a;
489  }
490  const double x = fabs(creal(z));
491  const double y = cimag(z), ya = fabs(y);
492 
493  cmplx ret = 0.; // return value
494 
495  double sum1 = 0, sum2 = 0, sum3 = 0, sum4 = 0, sum5 = 0;
496 
497  if ( ya > 7 || (x > 6
498  && (ya > 0.1 || (x > 8 && ya > 1e-10) || x > 28)) ) {
499 
500  // Poppe & Wijers suggest using a number of terms
501  // nu = 3 + 1442 / (26*rho + 77)
502  // where rho = sqrt((x/x0)^2 + (y/y0)^2) where x0=6.3, y0=4.4.
503  // (They only use this expansion for rho >= 1, but rho a little less
504  // than 1 seems okay too.)
505  // Instead, I did my own fit to a slightly different function
506  // that avoids the hypotenuse calculation, using NLopt to minimize
507  // the sum of the squares of the errors in nu with the constraint
508  // that the estimated nu be >= minimum nu to attain machine precision.
509  // I also separate the regions where nu == 2 and nu == 1.
510  const double ispi = 0.56418958354775628694807945156; // 1 / sqrt(pi)
511  double xs = y < 0 ? -creal(z) : creal(z); // compute for -z if y < 0
512  if ( x + ya > 4000 ) { // nu <= 2
513  if ( x + ya > 1e7 ) { // nu == 1, w(z) = i/sqrt(pi) / z
514  // scale to avoid overflow
515  if ( x > ya ) {
516  double yax = ya / xs;
517  double denom = ispi / (xs + yax*ya);
518  ret = C(denom*yax, denom);
519  } else if ( utility::is_inf(ya) ) {
520  return ((utility::is_nan(x) || y < 0) ? C(NaN,NaN) : C(0,0));
521  } else {
522  double xya = xs / ya;
523  double denom = ispi / (xya*xs + ya);
524  ret = C(denom, denom*xya);
525  }
526  } else { // nu == 2, w(z) = i/sqrt(pi) * z / (z*z - 0.5)
527  double dr = xs*xs - ya*ya - 0.5, di = 2*xs*ya;
528  double denom = ispi / (dr*dr + di*di);
529  ret = C(denom * (xs*di-ya*dr), denom * (xs*dr+ya*di));
530  }
531  } else { // compute nu(z) estimate and do general continued fraction
532  const double c0=3.9, c1=11.398, c2=0.08254, c3=0.1421, c4=0.2023; // fit
533  double nu = floor(c0 + c1 / (c2*x + c3*ya + c4));
534  double wr = xs, wi = ya;
535  for ( nu = 0.5 * (nu - 1); nu > 0.4; nu -= 0.5 ) {
536  // w <- z - nu/w:
537  double denom = nu / (wr*wr + wi*wi);
538  wr = xs - wr * denom;
539  wi = ya + wi * denom;
540  }
541  { // w(z) = i/sqrt(pi) / w:
542  double denom = ispi / (wr*wr + wi*wi);
543  ret = C(denom*wi, denom*wr);
544  }
545  }
546  if ( y < 0 ) {
547  // use w(z) = 2.0*exp(-z*z) - w(-z),
548  // but be careful of overflow in exp(-z*z)
549  // = exp(-(xs*xs-ya*ya) -2*i*xs*ya)
550  return 2.0*cexp(C((ya-xs)*(xs+ya), 2*xs*y)) - ret;
551  } else {
552  return ret;
553  }
554  } else if ( x < 10 ) {
555  // Note: The test that seems to be suggested in the paper is x <
556  // sqrt(-log(DBL_MIN)), about 26.6, since otherwise exp(-x^2)
557  // underflows to zero and sum1,sum2,sum4 are zero. However, long
558  // before this occurs, the sum1,sum2,sum4 contributions are
559  // negligible in double precision; I find that this happens for x >
560  // about 6, for all y. On the other hand, I find that the case
561  // where we compute all of the sums is faster (at least with the
562  // precomputed expa2n2 table) until about x=10. Furthermore, if we
563  // try to compute all of the sums for x > 20, I find that we
564  // sometimes run into numerical problems because underflow/overflow
565  // problems start to appear in the various coefficients of the sums,
566  // below. Therefore, we use x < 10 here.
567  double prod2ax = 1, prodm2ax = 1;
568  double expx2;
569 
570  if ( utility::is_nan(y) ) {
571  return C(y,y);
572  }
573 
574  // Somewhat ugly copy-and-paste duplication here, but I see significant
575  // speedups from using the special-case code with the precomputed
576  // exponential, and the x < 5e-4 special case is needed for accuracy.
577 
578  if ( relerr == DBL_EPSILON ) { // use precomputed exp(-a2*(n*n)) table
579  if ( x < 5e-4 ) { // compute sum4 and sum5 together as sum5-sum4
580  const double x2 = x*x;
581  expx2 = 1 - x2 * (1 - 0.5*x2); // exp(-x*x) via Taylor
582  // compute exp(2*a*x) and exp(-2*a*x) via Taylor, to double precision
583  const double ax2 = 1.036642960860171859744*x; // 2*a*x
584  const double exp2ax =
585  1 + ax2 * (1 + ax2 * (0.5 + 0.166666666666666666667*ax2));
586  const double expm2ax =
587  1 - ax2 * (1 - ax2 * (0.5 - 0.166666666666666666667*ax2));
588  for ( int n = 1; n <= int(expa2n2_length); ++n ) {
589  const double coef = expa2n2[n-1] * expx2 / (a2*(n*n) + y*y);
590  prod2ax *= exp2ax;
591  prodm2ax *= expm2ax;
592  sum1 += coef;
593  sum2 += coef * prodm2ax;
594  sum3 += coef * prod2ax;
595 
596  // really = sum5 - sum4
597  sum5 += coef * (2*a) * n * sinh_taylor((2*a)*n*x);
598 
599  // test convergence via sum3
600  if ( coef * prod2ax < relerr * sum3 ) break;
601  }
602  } else { // x > 5e-4, compute sum4 and sum5 separately
603  expx2 = exp(-x*x);
604  const double exp2ax = exp((2*a)*x), expm2ax = 1 / exp2ax;
605  for ( int n = 1; n <= int(expa2n2_length); ++n ) {
606  const double coef = expa2n2[n-1] * expx2 / (a2*(n*n) + y*y);
607  prod2ax *= exp2ax;
608  prodm2ax *= expm2ax;
609  sum1 += coef;
610  sum2 += coef * prodm2ax;
611  sum4 += (coef * prodm2ax) * (a*n);
612  sum3 += coef * prod2ax;
613  sum5 += (coef * prod2ax) * (a*n);
614  // test convergence via sum5, since this sum has the slowest decay
615  if ( (coef * prod2ax) * (a*n) < relerr * sum5 ) break;
616  }
617  }
618  } else { // relerr != DBL_EPSILON, compute exp(-a2*(n*n)) on the fly
619  const double exp2ax = exp((2*a)*x), expm2ax = 1 / exp2ax;
620  if ( x < 5e-4 ) { // compute sum4 and sum5 together as sum5-sum4
621  const double x2 = x*x;
622  expx2 = 1 - x2 * (1 - 0.5*x2); // exp(-x*x) via Taylor
623  for ( int n = 1; 1; ++n ) {
624  const double coef = exp(-a2*(n*n)) * expx2 / (a2*(n*n) + y*y);
625  prod2ax *= exp2ax;
626  prodm2ax *= expm2ax;
627  sum1 += coef;
628  sum2 += coef * prodm2ax;
629  sum3 += coef * prod2ax;
630 
631  // really = sum5 - sum4
632  sum5 += coef * (2*a) * n * sinh_taylor((2*a)*n*x);
633 
634  // test convergence via sum3
635  if ( coef * prod2ax < relerr * sum3 ) break;
636  }
637  } else { // x > 5e-4, compute sum4 and sum5 separately
638  expx2 = exp(-x*x);
639  for ( int n = 1; 1; ++n ) {
640  const double coef = exp(-a2*(n*n)) * expx2 / (a2*(n*n) + y*y);
641  prod2ax *= exp2ax;
642  prodm2ax *= expm2ax;
643  sum1 += coef;
644  sum2 += coef * prodm2ax;
645  sum4 += (coef * prodm2ax) * (a*n);
646  sum3 += coef * prod2ax;
647  sum5 += (coef * prod2ax) * (a*n);
648  // test convergence via sum5, since this sum has the slowest decay
649  if ( (coef * prod2ax) * (a*n) < relerr * sum5 ) break;
650  }
651  }
652  }
653  const double expx2erfcxy = // avoid spurious overflow for large negative y
654  y > -6 // for y < -6, erfcx(y) = 2*exp(y*y) to double precision
655  ? expx2*errfcx(y) : 2*exp(y*y-x*x);
656  if ( y > 5 ) { // imaginary terms cancel
657  const double sinxy = sin(x*y);
658  ret = (expx2erfcxy - c*y*sum1) * cos(2*x*y)
659  + (c*x*expx2) * sinxy * sinc(x*y, sinxy);
660  } else {
661  double xs = creal(z);
662  const double sinxy = sin(xs*y);
663  const double sin2xy = sin(2*xs*y), cos2xy = cos(2*xs*y);
664  const double coef1 = expx2erfcxy - c*y*sum1;
665  const double coef2 = c*xs*expx2;
666  ret = C(coef1 * cos2xy + coef2 * sinxy * sinc(xs*y, sinxy),
667  coef2 * sinc(2*xs*y, sin2xy) - coef1 * sin2xy);
668  }
669  } else { // x large: only sum3 & sum5 contribute (see above note)
670  if ( utility::is_nan(x) ) {
671  return C(x,x);
672  }
673  if ( utility::is_nan(y) ) {
674  return C(y,y);
675  }
676 
677  ret = exp(-x*x); // |y| < 1e-10, so we only need exp(-x*x) term
678 
679  // (round instead of ceil as in original paper; note that x/a > 1 here)
680  double n0 = floor(x/a + 0.5); // sum in both directions, starting at n0
681  double dx = a*n0 - x;
682  sum3 = exp(-dx*dx) / (a2*(n0*n0) + y*y);
683  sum5 = a*n0 * sum3;
684  double exp1 = exp(4*a*dx), exp1dn = 1;
685  int dn;
686  for ( dn = 1; n0 - dn > 0; ++dn ) { // loop over n0-dn and n0+dn terms
687  double np = n0 + dn, nm = n0 - dn;
688  double tp = exp(-sqr(a*dn+dx));
689  double tm = tp * (exp1dn *= exp1); // trick to get tm from tp
690  tp /= (a2*(np*np) + y*y);
691  tm /= (a2*(nm*nm) + y*y);
692  sum3 += tp + tm;
693  sum5 += a * (np * tp + nm * tm);
694  if ( a * (np * tp + nm * tm) < relerr * sum5 ) goto finish;
695  }
696  while ( 1 ) { // loop over n0+dn terms only (since n0-dn <= 0)
697  double np = n0 + dn++;
698  double tp = exp(-sqr(a*dn+dx)) / (a2*(np*np) + y*y);
699  sum3 += tp;
700  sum5 += a * np * tp;
701  if ( a * np * tp < relerr * sum5 ) goto finish;
702  }
703  }
704  finish:
705  return ret + C((0.5*c)*y*(sum2+sum3),
706  (0.5*c)*utility::copysign(sum5-sum4, creal(z)));
707 }
708 
709 /////////////////////////////////////////////////////////////////////////
710 
711 // erfcx(x) = exp(x^2) erfc(x) function, for real x, written by
712 // Steven G. Johnson, October 2012.
713 //
714 // This function combines a few different ideas.
715 //
716 // First, for x > 50, it uses a continued-fraction expansion (same as
717 // for the Faddeeva function, but with algebraic simplifications for z=i*x).
718 //
719 // Second, for 0 <= x <= 50, it uses Chebyshev polynomial approximations,
720 // but with two twists:
721 //
722 // a) It maps x to y = 4 / (4+x) in [0,1]. This simple transformation,
723 // inspired by a similar transformation in the octave-forge/specfun
724 // erfcx by Soren Hauberg, results in much faster Chebyshev convergence
725 // than other simple transformations I have examined.
726 //
727 // b) Instead of using a single Chebyshev polynomial for the entire
728 // [0,1] y interval, we break the interval up into 100 equal
729 // subintervals, with a switch/lookup table, and use much lower
730 // degree Chebyshev polynomials in each subinterval. This greatly
731 // improves performance in my tests.
732 //
733 // For x < 0, we use the relationship erfcx(-x) = 2 exp(x^2) - erfc(x),
734 // with the usual checks for overflow etcetera.
735 //
736 // Performance-wise, it seems to be substantially faster than either
737 // the SLATEC DERFC function [or an erfcx function derived therefrom]
738 // or Cody's CALERF function (from netlib.org/specfun), while
739 // retaining near machine precision in accuracy.
740 
741 // Given y100=100*y, where y = 4/(4+x) for x >= 0, compute erfc(x).
742 // Uses a look-up table of 100 different Chebyshev polynomials
743 // for y intervals [0,0.01], [0.01,0.02], ...., [0.99,1], generated
744 // with the help of Maple and a little shell script. This allows
745 // the Chebyshev polynomials to be of significantly lower degree (about 1/4)
746 // compared to fitting the whole [0,1] interval with a single polynomial.
747 static double erfcx_y100(double y100) {
748  switch ((int) y100) {
749  case 0 : {
750  double t = 2*y100 - 1;
751  return 0.70878032454106438663e-3 + (0.71234091047026302958e-3 + (0.35779077297597742384e-5 + (0.17403143962587937815e-7 + (0.81710660047307788845e-10 + (0.36885022360434957634e-12 + 0.15917038551111111111e-14 * t) * t) * t) * t) * t) * t;
752  }
753  case 1 : {
754  double t = 2*y100 - 3;
755  return 0.21479143208285144230e-2 + (0.72686402367379996033e-3 + (0.36843175430938995552e-5 + (0.18071841272149201685e-7 + (0.85496449296040325555e-10 + (0.38852037518534291510e-12 + 0.16868473576888888889e-14 * t) * t) * t) * t) * t) * t;
756  }
757  case 2 : {
758  double t = 2*y100 - 5;
759  return 0.36165255935630175090e-2 + (0.74182092323555510862e-3 + (0.37948319957528242260e-5 + (0.18771627021793087350e-7 + (0.89484715122415089123e-10 + (0.40935858517772440862e-12 + 0.17872061464888888889e-14 * t) * t) * t) * t) * t) * t;
760  }
761  case 3 : {
762  double t = 2*y100 - 7;
763  return 0.51154983860031979264e-2 + (0.75722840734791660540e-3 + (0.39096425726735703941e-5 + (0.19504168704300468210e-7 + (0.93687503063178993915e-10 + (0.43143925959079664747e-12 + 0.18939926435555555556e-14 * t) * t) * t) * t) * t) * t;
764  }
765  case 4 : {
766  double t = 2*y100 - 9;
767  return 0.66457513172673049824e-2 + (0.77310406054447454920e-3 + (0.40289510589399439385e-5 + (0.20271233238288381092e-7 + (0.98117631321709100264e-10 + (0.45484207406017752971e-12 + 0.20076352213333333333e-14 * t) * t) * t) * t) * t) * t;
768  }
769  case 5 : {
770  double t = 2*y100 - 11;
771  return 0.82082389970241207883e-2 + (0.78946629611881710721e-3 + (0.41529701552622656574e-5 + (0.21074693344544655714e-7 + (0.10278874108587317989e-9 + (0.47965201390613339638e-12 + 0.21285907413333333333e-14 * t) * t) * t) * t) * t) * t;
772  }
773  case 6 : {
774  double t = 2*y100 - 13;
775  return 0.98039537275352193165e-2 + (0.80633440108342840956e-3 + (0.42819241329736982942e-5 + (0.21916534346907168612e-7 + (0.10771535136565470914e-9 + (0.50595972623692822410e-12 + 0.22573462684444444444e-14 * t) * t) * t) * t) * t) * t;
776  }
777  case 7 : {
778  double t = 2*y100 - 15;
779  return 0.11433927298290302370e-1 + (0.82372858383196561209e-3 + (0.44160495311765438816e-5 + (0.22798861426211986056e-7 + (0.11291291745879239736e-9 + (0.53386189365816880454e-12 + 0.23944209546666666667e-14 * t) * t) * t) * t) * t) * t;
780  }
781  case 8 : {
782  double t = 2*y100 - 17;
783  return 0.13099232878814653979e-1 + (0.84167002467906968214e-3 + (0.45555958988457506002e-5 + (0.23723907357214175198e-7 + (0.11839789326602695603e-9 + (0.56346163067550237877e-12 + 0.25403679644444444444e-14 * t) * t) * t) * t) * t) * t;
784  }
785  case 9 : {
786  double t = 2*y100 - 19;
787  return 0.14800987015587535621e-1 + (0.86018092946345943214e-3 + (0.47008265848816866105e-5 + (0.24694040760197315333e-7 + (0.12418779768752299093e-9 + (0.59486890370320261949e-12 + 0.26957764568888888889e-14 * t) * t) * t) * t) * t) * t;
788  }
789  case 10 : {
790  double t = 2*y100 - 21;
791  return 0.16540351739394069380e-1 + (0.87928458641241463952e-3 + (0.48520195793001753903e-5 + (0.25711774900881709176e-7 + (0.13030128534230822419e-9 + (0.62820097586874779402e-12 + 0.28612737351111111111e-14 * t) * t) * t) * t) * t) * t;
792  }
793  case 11 : {
794  double t = 2*y100 - 23;
795  return 0.18318536789842392647e-1 + (0.89900542647891721692e-3 + (0.50094684089553365810e-5 + (0.26779777074218070482e-7 + (0.13675822186304615566e-9 + (0.66358287745352705725e-12 + 0.30375273884444444444e-14 * t) * t) * t) * t) * t) * t;
796  }
797  case 12 : {
798  double t = 2*y100 - 25;
799  return 0.20136801964214276775e-1 + (0.91936908737673676012e-3 + (0.51734830914104276820e-5 + (0.27900878609710432673e-7 + (0.14357976402809042257e-9 + (0.70114790311043728387e-12 + 0.32252476000000000000e-14 * t) * t) * t) * t) * t) * t;
800  }
801  case 13 : {
802  double t = 2*y100 - 27;
803  return 0.21996459598282740954e-1 + (0.94040248155366777784e-3 + (0.53443911508041164739e-5 + (0.29078085538049374673e-7 + (0.15078844500329731137e-9 + (0.74103813647499204269e-12 + 0.34251892320000000000e-14 * t) * t) * t) * t) * t) * t;
804  }
805  case 14 : {
806  double t = 2*y100 - 29;
807  return 0.23898877187226319502e-1 + (0.96213386835900177540e-3 + (0.55225386998049012752e-5 + (0.30314589961047687059e-7 + (0.15840826497296335264e-9 + (0.78340500472414454395e-12 + 0.36381553564444444445e-14 * t) * t) * t) * t) * t) * t;
808  }
809  case 15 : {
810  double t = 2*y100 - 31;
811  return 0.25845480155298518485e-1 + (0.98459293067820123389e-3 + (0.57082915920051843672e-5 + (0.31613782169164830118e-7 + (0.16646478745529630813e-9 + (0.82840985928785407942e-12 + 0.38649975768888888890e-14 * t) * t) * t) * t) * t) * t;
812  }
813  case 16 : {
814  double t = 2*y100 - 33;
815  return 0.27837754783474696598e-1 + (0.10078108563256892757e-2 + (0.59020366493792212221e-5 + (0.32979263553246520417e-7 + (0.17498524159268458073e-9 + (0.87622459124842525110e-12 + 0.41066206488888888890e-14 * t) * t) * t) * t) * t) * t;
816  }
817  case 17 : {
818  double t = 2*y100 - 35;
819  return 0.29877251304899307550e-1 + (0.10318204245057349310e-2 + (0.61041829697162055093e-5 + (0.34414860359542720579e-7 + (0.18399863072934089607e-9 + (0.92703227366365046533e-12 + 0.43639844053333333334e-14 * t) * t) * t) * t) * t) * t;
820  }
821  case 18 : {
822  double t = 2*y100 - 37;
823  return 0.31965587178596443475e-1 + (0.10566560976716574401e-2 + (0.63151633192414586770e-5 + (0.35924638339521924242e-7 + (0.19353584758781174038e-9 + (0.98102783859889264382e-12 + 0.46381060817777777779e-14 * t) * t) * t) * t) * t) * t;
824  }
825  case 19 : {
826  double t = 2*y100 - 39;
827  return 0.34104450552588334840e-1 + (0.10823541191350532574e-2 + (0.65354356159553934436e-5 + (0.37512918348533521149e-7 + (0.20362979635817883229e-9 + (0.10384187833037282363e-11 + 0.49300625262222222221e-14 * t) * t) * t) * t) * t) * t;
828  }
829  case 20 : {
830  double t = 2*y100 - 41;
831  return 0.36295603928292425716e-1 + (0.11089526167995268200e-2 + (0.67654845095518363577e-5 + (0.39184292949913591646e-7 + (0.21431552202133775150e-9 + (0.10994259106646731797e-11 + 0.52409949102222222221e-14 * t) * t) * t) * t) * t) * t;
832  }
833  case 21 : {
834  double t = 2*y100 - 43;
835  return 0.38540888038840509795e-1 + (0.11364917134175420009e-2 + (0.70058230641246312003e-5 + (0.40943644083718586939e-7 + (0.22563034723692881631e-9 + (0.11642841011361992885e-11 + 0.55721092871111111110e-14 * t) * t) * t) * t) * t) * t;
836  }
837  case 22 : {
838  double t = 2*y100 - 45;
839  return 0.40842225954785960651e-1 + (0.11650136437945673891e-2 + (0.72569945502343006619e-5 + (0.42796161861855042273e-7 + (0.23761401711005024162e-9 + (0.12332431172381557035e-11 + 0.59246802364444444445e-14 * t) * t) * t) * t) * t) * t;
840  }
841  case 23 : {
842  double t = 2*y100 - 47;
843  return 0.43201627431540222422e-1 + (0.11945628793917272199e-2 + (0.75195743532849206263e-5 + (0.44747364553960993492e-7 + (0.25030885216472953674e-9 + (0.13065684400300476484e-11 + 0.63000532853333333334e-14 * t) * t) * t) * t) * t) * t;
844  }
845  case 24 : {
846  double t = 2*y100 - 49;
847  return 0.45621193513810471438e-1 + (0.12251862608067529503e-2 + (0.77941720055551920319e-5 + (0.46803119830954460212e-7 + (0.26375990983978426273e-9 + (0.13845421370977119765e-11 + 0.66996477404444444445e-14 * t) * t) * t) * t) * t) * t;
848  }
849  case 25 : {
850  double t = 2*y100 - 51;
851  return 0.48103121413299865517e-1 + (0.12569331386432195113e-2 + (0.80814333496367673980e-5 + (0.48969667335682018324e-7 + (0.27801515481905748484e-9 + (0.14674637611609884208e-11 + 0.71249589351111111110e-14 * t) * t) * t) * t) * t) * t;
852  }
853  case 26 : {
854  double t = 2*y100 - 53;
855  return 0.50649709676983338501e-1 + (0.12898555233099055810e-2 + (0.83820428414568799654e-5 + (0.51253642652551838659e-7 + (0.29312563849675507232e-9 + (0.15556512782814827846e-11 + 0.75775607822222222221e-14 * t) * t) * t) * t) * t) * t;
856  }
857  case 27 : {
858  double t = 2*y100 - 55;
859  return 0.53263363664388864181e-1 + (0.13240082443256975769e-2 + (0.86967260015007658418e-5 + (0.53662102750396795566e-7 + (0.30914568786634796807e-9 + (0.16494420240828493176e-11 + 0.80591079644444444445e-14 * t) * t) * t) * t) * t) * t;
860  }
861  case 28 : {
862  double t = 2*y100 - 57;
863  return 0.55946601353500013794e-1 + (0.13594491197408190706e-2 + (0.90262520233016380987e-5 + (0.56202552975056695376e-7 + (0.32613310410503135996e-9 + (0.17491936862246367398e-11 + 0.85713381688888888890e-14 * t) * t) * t) * t) * t) * t;
864  }
865  case 29 : {
866  double t = 2*y100 - 59;
867  return 0.58702059496154081813e-1 + (0.13962391363223647892e-2 + (0.93714365487312784270e-5 + (0.58882975670265286526e-7 + (0.34414937110591753387e-9 + (0.18552853109751857859e-11 + 0.91160736711111111110e-14 * t) * t) * t) * t) * t) * t;
868  }
869  case 30 : {
870  double t = 2*y100 - 61;
871  return 0.61532500145144778048e-1 + (0.14344426411912015247e-2 + (0.97331446201016809696e-5 + (0.61711860507347175097e-7 + (0.36325987418295300221e-9 + (0.19681183310134518232e-11 + 0.96952238400000000000e-14 * t) * t) * t) * t) * t) * t;
872  }
873  case 31 : {
874  double t = 2*y100 - 63;
875  return 0.64440817576653297993e-1 + (0.14741275456383131151e-2 + (0.10112293819576437838e-4 + (0.64698236605933246196e-7 + (0.38353412915303665586e-9 + (0.20881176114385120186e-11 + 0.10310784480000000000e-13 * t) * t) * t) * t) * t) * t;
876  }
877  case 32 : {
878  double t = 2*y100 - 65;
879  return 0.67430045633130393282e-1 + (0.15153655418916540370e-2 + (0.10509857606888328667e-4 + (0.67851706529363332855e-7 + (0.40504602194811140006e-9 + (0.22157325110542534469e-11 + 0.10964842115555555556e-13 * t) * t) * t) * t) * t) * t;
880  }
881  case 33 : {
882  double t = 2*y100 - 67;
883  return 0.70503365513338850709e-1 + (0.15582323336495709827e-2 + (0.10926868866865231089e-4 + (0.71182482239613507542e-7 + (0.42787405890153386710e-9 + (0.23514379522274416437e-11 + 0.11659571751111111111e-13 * t) * t) * t) * t) * t) * t;
884  }
885  case 34 : {
886  double t = 2*y100 - 69;
887  return 0.73664114037944596353e-1 + (0.16028078812438820413e-2 + (0.11364423678778207991e-4 + (0.74701423097423182009e-7 + (0.45210162777476488324e-9 + (0.24957355004088569134e-11 + 0.12397238257777777778e-13 * t) * t) * t) * t) * t) * t;
888  }
889  case 35 : {
890  double t = 2*y100 - 71;
891  return 0.76915792420819562379e-1 + (0.16491766623447889354e-2 + (0.11823685320041302169e-4 + (0.78420075993781544386e-7 + (0.47781726956916478925e-9 + (0.26491544403815724749e-11 + 0.13180196462222222222e-13 * t) * t) * t) * t) * t) * t;
892  }
893  case 36 : {
894  double t = 2*y100 - 73;
895  return 0.80262075578094612819e-1 + (0.16974279491709504117e-2 + (0.12305888517309891674e-4 + (0.82350717698979042290e-7 + (0.50511496109857113929e-9 + (0.28122528497626897696e-11 + 0.14010889635555555556e-13 * t) * t) * t) * t) * t) * t;
896  }
897  case 37 : {
898  double t = 2*y100 - 75;
899  return 0.83706822008980357446e-1 + (0.17476561032212656962e-2 + (0.12812343958540763368e-4 + (0.86506399515036435592e-7 + (0.53409440823869467453e-9 + (0.29856186620887555043e-11 + 0.14891851591111111111e-13 * t) * t) * t) * t) * t) * t;
900  }
901  case 38 : {
902  double t = 2*y100 - 77;
903  return 0.87254084284461718231e-1 + (0.17999608886001962327e-2 + (0.13344443080089492218e-4 + (0.90900994316429008631e-7 + (0.56486134972616465316e-9 + (0.31698707080033956934e-11 + 0.15825697795555555556e-13 * t) * t) * t) * t) * t) * t;
904  }
905  case 39 : {
906  double t = 2*y100 - 79;
907  return 0.90908120182172748487e-1 + (0.18544478050657699758e-2 + (0.13903663143426120077e-4 + (0.95549246062549906177e-7 + (0.59752787125242054315e-9 + (0.33656597366099099413e-11 + 0.16815130613333333333e-13 * t) * t) * t) * t) * t) * t;
908  }
909  case 40 : {
910  double t = 2*y100 - 81;
911  return 0.94673404508075481121e-1 + (0.19112284419887303347e-2 + (0.14491572616545004930e-4 + (0.10046682186333613697e-6 + (0.63221272959791000515e-9 + (0.35736693975589130818e-11 + 0.17862931591111111111e-13 * t) * t) * t) * t) * t) * t;
912  }
913  case 41 : {
914  double t = 2*y100 - 83;
915  return 0.98554641648004456555e-1 + (0.19704208544725622126e-2 + (0.15109836875625443935e-4 + (0.10567036667675984067e-6 + (0.66904168640019354565e-9 + (0.37946171850824333014e-11 + 0.18971959040000000000e-13 * t) * t) * t) * t) * t) * t;
916  }
917  case 42 : {
918  double t = 2*y100 - 85;
919  return 0.10255677889470089531e0 + (0.20321499629472857418e-2 + (0.15760224242962179564e-4 + (0.11117756071353507391e-6 + (0.70814785110097658502e-9 + (0.40292553276632563925e-11 + 0.20145143075555555556e-13 * t) * t) * t) * t) * t) * t;
920  }
921  case 43 : {
922  double t = 2*y100 - 87;
923  return 0.10668502059865093318e0 + (0.20965479776148731610e-2 + (0.16444612377624983565e-4 + (0.11700717962026152749e-6 + (0.74967203250938418991e-9 + (0.42783716186085922176e-11 + 0.21385479360000000000e-13 * t) * t) * t) * t) * t) * t;
924  }
925  case 44 : {
926  double t = 2*y100 - 89;
927  return 0.11094484319386444474e0 + (0.21637548491908170841e-2 + (0.17164995035719657111e-4 + (0.12317915750735938089e-6 + (0.79376309831499633734e-9 + (0.45427901763106353914e-11 + 0.22696025653333333333e-13 * t) * t) * t) * t) * t) * t;
928  }
929  case 45 : {
930  double t = 2*y100 - 91;
931  return 0.11534201115268804714e0 + (0.22339187474546420375e-2 + (0.17923489217504226813e-4 + (0.12971465288245997681e-6 + (0.84057834180389073587e-9 + (0.48233721206418027227e-11 + 0.24079890062222222222e-13 * t) * t) * t) * t) * t) * t;
932  }
933  case 46 : {
934  double t = 2*y100 - 93;
935  return 0.11988259392684094740e0 + (0.23071965691918689601e-2 + (0.18722342718958935446e-4 + (0.13663611754337957520e-6 + (0.89028385488493287005e-9 + (0.51210161569225846701e-11 + 0.25540227111111111111e-13 * t) * t) * t) * t) * t) * t;
936  }
937  case 47 : {
938  double t = 2*y100 - 95;
939  return 0.12457298393509812907e0 + (0.23837544771809575380e-2 + (0.19563942105711612475e-4 + (0.14396736847739470782e-6 + (0.94305490646459247016e-9 + (0.54366590583134218096e-11 + 0.27080225920000000000e-13 * t) * t) * t) * t) * t) * t;
940  }
941  case 48 : {
942  double t = 2*y100 - 97;
943  return 0.12941991566142438816e0 + (0.24637684719508859484e-2 + (0.20450821127475879816e-4 + (0.15173366280523906622e-6 + (0.99907632506389027739e-9 + (0.57712760311351625221e-11 + 0.28703099555555555556e-13 * t) * t) * t) * t) * t) * t;
944  }
945  case 49 : {
946  double t = 2*y100 - 99;
947  return 0.13443048593088696613e0 + (0.25474249981080823877e-2 + (0.21385669591362915223e-4 + (0.15996177579900443030e-6 + (0.10585428844575134013e-8 + (0.61258809536787882989e-11 + 0.30412080142222222222e-13 * t) * t) * t) * t) * t) * t;
948  }
949  case 50 : {
950  double t = 2*y100 - 101;
951  return 0.13961217543434561353e0 + (0.26349215871051761416e-2 + (0.22371342712572567744e-4 + (0.16868008199296822247e-6 + (0.11216596910444996246e-8 + (0.65015264753090890662e-11 + 0.32210394506666666666e-13 * t) * t) * t) * t) * t) * t;
952  }
953  case 51 : {
954  double t = 2*y100 - 103;
955  return 0.14497287157673800690e0 + (0.27264675383982439814e-2 + (0.23410870961050950197e-4 + (0.17791863939526376477e-6 + (0.11886425714330958106e-8 + (0.68993039665054288034e-11 + 0.34101266222222222221e-13 * t) * t) * t) * t) * t) * t;
956  }
957  case 52 : {
958  double t = 2*y100 - 105;
959  return 0.15052089272774618151e0 + (0.28222846410136238008e-2 + (0.24507470422713397006e-4 + (0.18770927679626136909e-6 + (0.12597184587583370712e-8 + (0.73203433049229821618e-11 + 0.36087889048888888890e-13 * t) * t) * t) * t) * t) * t;
960  }
961  case 53 : {
962  double t = 2*y100 - 107;
963  return 0.15626501395774612325e0 + (0.29226079376196624949e-2 + (0.25664553693768450545e-4 + (0.19808568415654461964e-6 + (0.13351257759815557897e-8 + (0.77658124891046760667e-11 + 0.38173420035555555555e-13 * t) * t) * t) * t) * t) * t;
964  }
965  case 54 : {
966  double t = 2*y100 - 109;
967  return 0.16221449434620737567e0 + (0.30276865332726475672e-2 + (0.26885741326534564336e-4 + (0.20908350604346384143e-6 + (0.14151148144240728728e-8 + (0.82369170665974313027e-11 + 0.40360957457777777779e-13 * t) * t) * t) * t) * t) * t;
968  }
969  case 55 : {
970  double t = 2*y100 - 111;
971  return 0.16837910595412130659e0 + (0.31377844510793082301e-2 + (0.28174873844911175026e-4 + (0.22074043807045782387e-6 + (0.14999481055996090039e-8 + (0.87348993661930809254e-11 + 0.42653528977777777779e-13 * t) * t) * t) * t) * t) * t;
972  }
973  case 56 : {
974  double t = 2*y100 - 113;
975  return 0.17476916455659369953e0 + (0.32531815370903068316e-2 + (0.29536024347344364074e-4 + (0.23309632627767074202e-6 + (0.15899007843582444846e-8 + (0.92610375235427359475e-11 + 0.45054073102222222221e-13 * t) * t) * t) * t) * t) * t;
976  }
977  case 57 : {
978  double t = 2*y100 - 115;
979  return 0.18139556223643701364e0 + (0.33741744168096996041e-2 + (0.30973511714709500836e-4 + (0.24619326937592290996e-6 + (0.16852609412267750744e-8 + (0.98166442942854895573e-11 + 0.47565418097777777779e-13 * t) * t) * t) * t) * t) * t;
980  }
981  case 58 : {
982  double t = 2*y100 - 117;
983  return 0.18826980194443664549e0 + (0.35010775057740317997e-2 + (0.32491914440014267480e-4 + (0.26007572375886319028e-6 + (0.17863299617388376116e-8 + (0.10403065638343878679e-10 + 0.50190265831111111110e-13 * t) * t) * t) * t) * t) * t;
984  }
985  case 59 : {
986  double t = 2*y100 - 119;
987  return 0.19540403413693967350e0 + (0.36342240767211326315e-2 + (0.34096085096200907289e-4 + (0.27479061117017637474e-6 + (0.18934228504790032826e-8 + (0.11021679075323598664e-10 + 0.52931171733333333334e-13 * t) * t) * t) * t) * t) * t;
988  }
989  case 60 : {
990  double t = 2*y100 - 121;
991  return 0.20281109560651886959e0 + (0.37739673859323597060e-2 + (0.35791165457592409054e-4 + (0.29038742889416172404e-6 + (0.20068685374849001770e-8 + (0.11673891799578381999e-10 + 0.55790523093333333334e-13 * t) * t) * t) * t) * t) * t;
992  }
993  case 61 : {
994  double t = 2*y100 - 123;
995  return 0.21050455062669334978e0 + (0.39206818613925652425e-2 + (0.37582602289680101704e-4 + (0.30691836231886877385e-6 + (0.21270101645763677824e-8 + (0.12361138551062899455e-10 + 0.58770520160000000000e-13 * t) * t) * t) * t) * t) * t;
996  }
997  case 62 : {
998  double t = 2*y100 - 125;
999  return 0.21849873453703332479e0 + (0.40747643554689586041e-2 + (0.39476163820986711501e-4 + (0.32443839970139918836e-6 + (0.22542053491518680200e-8 + (0.13084879235290858490e-10 + 0.61873153262222222221e-13 * t) * t) * t) * t) * t) * t;
1000  }
1001  case 63 : {
1002  double t = 2*y100 - 127;
1003  return 0.22680879990043229327e0 + (0.42366354648628516935e-2 + (0.41477956909656896779e-4 + (0.34300544894502810002e-6 + (0.23888264229264067658e-8 + (0.13846596292818514601e-10 + 0.65100183751111111110e-13 * t) * t) * t) * t) * t) * t;
1004  }
1005  case 64 : {
1006  double t = 2*y100 - 129;
1007  return 0.23545076536988703937e0 + (0.44067409206365170888e-2 + (0.43594444916224700881e-4 + (0.36268045617760415178e-6 + (0.25312606430853202748e-8 + (0.14647791812837903061e-10 + 0.68453122631111111110e-13 * t) * t) * t) * t) * t) * t;
1008  }
1009  case 65 : {
1010  double t = 2*y100 - 131;
1011  return 0.24444156740777432838e0 + (0.45855530511605787178e-2 + (0.45832466292683085475e-4 + (0.38352752590033030472e-6 + (0.26819103733055603460e-8 + (0.15489984390884756993e-10 + 0.71933206364444444445e-13 * t) * t) * t) * t) * t) * t;
1012  }
1013  case 66 : {
1014  double t = 2*y100 - 133;
1015  return 0.25379911500634264643e0 + (0.47735723208650032167e-2 + (0.48199253896534185372e-4 + (0.40561404245564732314e-6 + (0.28411932320871165585e-8 + (0.16374705736458320149e-10 + 0.75541379822222222221e-13 * t) * t) * t) * t) * t) * t;
1016  }
1017  case 67 : {
1018  double t = 2*y100 - 135;
1019  return 0.26354234756393613032e0 + (0.49713289477083781266e-2 + (0.50702455036930367504e-4 + (0.42901079254268185722e-6 + (0.30095422058900481753e-8 + (0.17303497025347342498e-10 + 0.79278273368888888890e-13 * t) * t) * t) * t) * t) * t;
1020  }
1021  case 68 : {
1022  double t = 2*y100 - 137;
1023  return 0.27369129607732343398e0 + (0.51793846023052643767e-2 + (0.53350152258326602629e-4 + (0.45379208848865015485e-6 + (0.31874057245814381257e-8 + (0.18277905010245111046e-10 + 0.83144182364444444445e-13 * t) * t) * t) * t) * t) * t;
1024  }
1025  case 69 : {
1026  double t = 2*y100 - 139;
1027  return 0.28426714781640316172e0 + (0.53983341916695141966e-2 + (0.56150884865255810638e-4 + (0.48003589196494734238e-6 + (0.33752476967570796349e-8 + (0.19299477888083469086e-10 + 0.87139049137777777779e-13 * t) * t) * t) * t) * t) * t;
1028  }
1029  case 70 : {
1030  double t = 2*y100 - 141;
1031  return 0.29529231465348519920e0 + (0.56288077305420795663e-2 + (0.59113671189913307427e-4 + (0.50782393781744840482e-6 + (0.35735475025851713168e-8 + (0.20369760937017070382e-10 + 0.91262442613333333334e-13 * t) * t) * t) * t) * t) * t;
1032  }
1033  case 71 : {
1034  double t = 2*y100 - 143;
1035  return 0.30679050522528838613e0 + (0.58714723032745403331e-2 + (0.62248031602197686791e-4 + (0.53724185766200945789e-6 + (0.37827999418960232678e-8 + (0.21490291930444538307e-10 + 0.95513539182222222221e-13 * t) * t) * t) * t) * t) * t;
1036  }
1037  case 72 : {
1038  double t = 2*y100 - 145;
1039  return 0.31878680111173319425e0 + (0.61270341192339103514e-2 + (0.65564012259707640976e-4 + (0.56837930287837738996e-6 + (0.40035151353392378882e-8 + (0.22662596341239294792e-10 + 0.99891109760000000000e-13 * t) * t) * t) * t) * t) * t;
1040  }
1041  case 73 : {
1042  double t = 2*y100 - 147;
1043  return 0.33130773722152622027e0 + (0.63962406646798080903e-2 + (0.69072209592942396666e-4 + (0.60133006661885941812e-6 + (0.42362183765883466691e-8 + (0.23888182347073698382e-10 + 0.10439349811555555556e-12 * t) * t) * t) * t) * t) * t;
1044  }
1045  case 74 : {
1046  double t = 2*y100 - 149;
1047  return 0.34438138658041336523e0 + (0.66798829540414007258e-2 + (0.72783795518603561144e-4 + (0.63619220443228800680e-6 + (0.44814499336514453364e-8 + (0.25168535651285475274e-10 + 0.10901861383111111111e-12 * t) * t) * t) * t) * t) * t;
1048  }
1049  case 75 : {
1050  double t = 2*y100 - 151;
1051  return 0.35803744972380175583e0 + (0.69787978834882685031e-2 + (0.76710543371454822497e-4 + (0.67306815308917386747e-6 + (0.47397647975845228205e-8 + (0.26505114141143050509e-10 + 0.11376390933333333333e-12 * t) * t) * t) * t) * t) * t;
1052  }
1053  case 76 : {
1054  double t = 2*y100 - 153;
1055  return 0.37230734890119724188e0 + (0.72938706896461381003e-2 + (0.80864854542670714092e-4 + (0.71206484718062688779e-6 + (0.50117323769745883805e-8 + (0.27899342394100074165e-10 + 0.11862637614222222222e-12 * t) * t) * t) * t) * t) * t;
1056  }
1057  case 77 : {
1058  double t = 2*y100 - 155;
1059  return 0.38722432730555448223e0 + (0.76260375162549802745e-2 + (0.85259785810004603848e-4 + (0.75329383305171327677e-6 + (0.52979361368388119355e-8 + (0.29352606054164086709e-10 + 0.12360253370666666667e-12 * t) * t) * t) * t) * t) * t;
1060  }
1061  case 78 : {
1062  double t = 2*y100 - 157;
1063  return 0.40282355354616940667e0 + (0.79762880915029728079e-2 + (0.89909077342438246452e-4 + (0.79687137961956194579e-6 + (0.55989731807360403195e-8 + (0.30866246101464869050e-10 + 0.12868841946666666667e-12 * t) * t) * t) * t) * t) * t;
1064  }
1065  case 79 : {
1066  double t = 2*y100 - 159;
1067  return 0.41914223158913787649e0 + (0.83456685186950463538e-2 + (0.94827181359250161335e-4 + (0.84291858561783141014e-6 + (0.59154537751083485684e-8 + (0.32441553034347469291e-10 + 0.13387957943111111111e-12 * t) * t) * t) * t) * t) * t;
1068  }
1069  case 80 : {
1070  double t = 2*y100 - 161;
1071  return 0.43621971639463786896e0 + (0.87352841828289495773e-2 + (0.10002929142066799966e-3 + (0.89156148280219880024e-6 + (0.62480008150788597147e-8 + (0.34079760983458878910e-10 + 0.13917107176888888889e-12 * t) * t) * t) * t) * t) * t;
1072  }
1073  case 81 : {
1074  double t = 2*y100 - 163;
1075  return 0.45409763548534330981e0 + (0.91463027755548240654e-2 + (0.10553137232446167258e-3 + (0.94293113464638623798e-6 + (0.65972492312219959885e-8 + (0.35782041795476563662e-10 + 0.14455745872000000000e-12 * t) * t) * t) * t) * t) * t;
1076  }
1077  case 82 : {
1078  double t = 2*y100 - 165;
1079  return 0.47282001668512331468e0 + (0.95799574408860463394e-2 + (0.11135019058000067469e-3 + (0.99716373005509038080e-6 + (0.69638453369956970347e-8 + (0.37549499088161345850e-10 + 0.15003280712888888889e-12 * t) * t) * t) * t) * t) * t;
1080  }
1081  case 83 : {
1082  double t = 2*y100 - 167;
1083  return 0.49243342227179841649e0 + (0.10037550043909497071e-1 + (0.11750334542845234952e-3 + (0.10544006716188967172e-5 + (0.73484461168242224872e-8 + (0.39383162326435752965e-10 + 0.15559069118222222222e-12 * t) * t) * t) * t) * t) * t;
1084  }
1085  case 84 : {
1086  double t = 2*y100 - 169;
1087  return 0.51298708979209258326e0 + (0.10520454564612427224e-1 + (0.12400930037494996655e-3 + (0.11147886579371265246e-5 + (0.77517184550568711454e-8 + (0.41283980931872622611e-10 + 0.16122419680000000000e-12 * t) * t) * t) * t) * t) * t;
1088  }
1089  case 85 : {
1090  double t = 2*y100 - 171;
1091  return 0.53453307979101369843e0 + (0.11030120618800726938e-1 + (0.13088741519572269581e-3 + (0.11784797595374515432e-5 + (0.81743383063044825400e-8 + (0.43252818449517081051e-10 + 0.16692592640000000000e-12 * t) * t) * t) * t) * t) * t;
1092  }
1093  case 86 : {
1094  double t = 2*y100 - 173;
1095  return 0.55712643071169299478e0 + (0.11568077107929735233e-1 + (0.13815797838036651289e-3 + (0.12456314879260904558e-5 + (0.86169898078969313597e-8 + (0.45290446811539652525e-10 + 0.17268801084444444444e-12 * t) * t) * t) * t) * t) * t;
1096  }
1097  case 87 : {
1098  double t = 2*y100 - 175;
1099  return 0.58082532122519320968e0 + (0.12135935999503877077e-1 + (0.14584223996665838559e-3 + (0.13164068573095710742e-5 + (0.90803643355106020163e-8 + (0.47397540713124619155e-10 + 0.17850211608888888889e-12 * t) * t) * t) * t) * t) * t;
1100  }
1101  case 88 : {
1102  double t = 2*y100 - 177;
1103  return 0.60569124025293375554e0 + (0.12735396239525550361e-1 + (0.15396244472258863344e-3 + (0.13909744385382818253e-5 + (0.95651595032306228245e-8 + (0.49574672127669041550e-10 + 0.18435945564444444444e-12 * t) * t) * t) * t) * t) * t;
1104  }
1105  case 89 : {
1106  double t = 2*y100 - 179;
1107  return 0.63178916494715716894e0 + (0.13368247798287030927e-1 + (0.16254186562762076141e-3 + (0.14695084048334056083e-5 + (0.10072078109604152350e-7 + (0.51822304995680707483e-10 + 0.19025081422222222222e-12 * t) * t) * t) * t) * t) * t;
1108  }
1109  case 90 : {
1110  double t = 2*y100 - 181;
1111  return 0.65918774689725319200e0 + (0.14036375850601992063e-1 + (0.17160483760259706354e-3 + (0.15521885688723188371e-5 + (0.10601827031535280590e-7 + (0.54140790105837520499e-10 + 0.19616655146666666667e-12 * t) * t) * t) * t) * t) * t;
1112  }
1113  case 91 : {
1114  double t = 2*y100 - 183;
1115  return 0.68795950683174433822e0 + (0.14741765091365869084e-1 + (0.18117679143520433835e-3 + (0.16392004108230585213e-5 + (0.11155116068018043001e-7 + (0.56530360194925690374e-10 + 0.20209663662222222222e-12 * t) * t) * t) * t) * t) * t;
1116  }
1117  case 92 : {
1118  double t = 2*y100 - 185;
1119  return 0.71818103808729967036e0 + (0.15486504187117112279e-1 + (0.19128428784550923217e-3 + (0.17307350969359975848e-5 + (0.11732656736113607751e-7 + (0.58991125287563833603e-10 + 0.20803065333333333333e-12 * t) * t) * t) * t) * t) * t;
1120  }
1121  case 93 : {
1122  double t = 2*y100 - 187;
1123  return 0.74993321911726254661e0 + (0.16272790364044783382e-1 + (0.20195505163377912645e-3 + (0.18269894883203346953e-5 + (0.12335161021630225535e-7 + (0.61523068312169087227e-10 + 0.21395783431111111111e-12 * t) * t) * t) * t) * t) * t;
1124  }
1125  case 94 : {
1126  double t = 2*y100 - 189;
1127  return 0.78330143531283492729e0 + (0.17102934132652429240e-1 + (0.21321800585063327041e-3 + (0.19281661395543913713e-5 + (0.12963340087354341574e-7 + (0.64126040998066348872e-10 + 0.21986708942222222222e-12 * t) * t) * t) * t) * t) * t;
1128  }
1129  case 95 : {
1130  double t = 2*y100 - 191;
1131  return 0.81837581041023811832e0 + (0.17979364149044223802e-1 + (0.22510330592753129006e-3 + (0.20344732868018175389e-5 + (0.13617902941839949718e-7 + (0.66799760083972474642e-10 + 0.22574701262222222222e-12 * t) * t) * t) * t) * t) * t;
1132  }
1133  case 96 : {
1134  double t = 2*y100 - 193;
1135  return 0.85525144775685126237e0 + (0.18904632212547561026e-1 + (0.23764237370371255638e-3 + (0.21461248251306387979e-5 + (0.14299555071870523786e-7 + (0.69543803864694171934e-10 + 0.23158593688888888889e-12 * t) * t) * t) * t) * t) * t;
1136  }
1137  case 97 : {
1138  double t = 2*y100 - 195;
1139  return 0.89402868170849933734e0 + (0.19881418399127202569e-1 + (0.25086793128395995798e-3 + (0.22633402747585233180e-5 + (0.15008997042116532283e-7 + (0.72357609075043941261e-10 + 0.23737194737777777778e-12 * t) * t) * t) * t) * t) * t;
1140  }
1141  case 98 : {
1142  double t = 2*y100 - 197;
1143  return 0.93481333942870796363e0 + (0.20912536329780368893e-1 + (0.26481403465998477969e-3 + (0.23863447359754921676e-5 + (0.15746923065472184451e-7 + (0.75240468141720143653e-10 + 0.24309291271111111111e-12 * t) * t) * t) * t) * t) * t;
1144  }
1145  case 99 : {
1146  double t = 2*y100 - 199;
1147  return 0.97771701335885035464e0 + (0.22000938572830479551e-1 + (0.27951610702682383001e-3 + (0.25153688325245314530e-5 + (0.16514019547822821453e-7 + (0.78191526829368231251e-10 + 0.24873652355555555556e-12 * t) * t) * t) * t) * t) * t;
1148  }
1149  }
1150  // we only get here if y = 1, i.e. |x| < 4*eps, in which case
1151  // erfcx is within 1e-15 of 1..
1152  return 1.0;
1153 }
1154 
1155 double errfcx(double x) {
1156  if ( x >= 0 ) {
1157  if ( x > 50 ) { // continued-fraction expansion is faster
1158  const double ispi = 0.56418958354775628694807945156; // 1 / sqrt(pi)
1159  if ( x > 5e7 ) { // 1-term expansion, important to avoid overflow
1160  return ispi / x;
1161  }
1162  // 5-term expansion (rely on compiler for CSE), simplified from: ispi / (x+0.5/(x+1/(x+1.5/(x+2/x))))
1163  return ispi*((x*x) * (x*x+4.5) + 2) / (x * ((x*x) * (x*x+5) + 3.75));
1164  }
1165  return erfcx_y100(400/(4+x));
1166  } else {
1167  return x < -26.7 ? HUGE_VAL : (x < -6.1 ? 2*exp(x*x)
1168  : 2*exp(x*x) - erfcx_y100(400/(4-x)));
1169  }
1170 }
1171 
1172 /////////////////////////////////////////////////////////////////////////
1173 // Compute a scaled Dawson integral
1174 // w_im(x) = 2*Dawson(x)/sqrt(pi)
1175 // equivalent to the imaginary part w(x) for real x.
1176 //
1177 // Uses methods similar to the erfcx calculation above: continued fractions
1178 // for large |x|, a lookup table of Chebyshev polynomials for smaller |x|,
1179 // and finally a Taylor expansion for |x|<0.01.
1180 //
1181 // Steven G. Johnson, October 2012.
1182 
1183 // Given y100=100*y, where y = 1/(1+x) for x >= 0, compute w_im(x).
1184 //
1185 // Uses a look-up table of 100 different Chebyshev polynomials
1186 // for y intervals [0,0.01], [0.01,0.02], ...., [0.99,1], generated
1187 // with the help of Maple and a little shell script. This allows
1188 // the Chebyshev polynomials to be of significantly lower degree (about 1/30)
1189 // compared to fitting the whole [0,1] interval with a single polynomial.
1190 static double w_im_y100(double y100, double x) {
1191  switch ((int) y100) {
1192  case 0 : {
1193  double t = 2*y100 - 1;
1194  return 0.28351593328822191546e-2 + (0.28494783221378400759e-2 + (0.14427470563276734183e-4 + (0.10939723080231588129e-6 + (0.92474307943275042045e-9 + (0.89128907666450075245e-11 + 0.92974121935111111110e-13 * t) * t) * t) * t) * t) * t;
1195  }
1196  case 1 : {
1197  double t = 2*y100 - 3;
1198  return 0.85927161243940350562e-2 + (0.29085312941641339862e-2 + (0.15106783707725582090e-4 + (0.11716709978531327367e-6 + (0.10197387816021040024e-8 + (0.10122678863073360769e-10 + 0.10917479678400000000e-12 * t) * t) * t) * t) * t) * t;
1199  }
1200  case 2 : {
1201  double t = 2*y100 - 5;
1202  return 0.14471159831187703054e-1 + (0.29703978970263836210e-2 + (0.15835096760173030976e-4 + (0.12574803383199211596e-6 + (0.11278672159518415848e-8 + (0.11547462300333495797e-10 + 0.12894535335111111111e-12 * t) * t) * t) * t) * t) * t;
1203  }
1204  case 3 : {
1205  double t = 2*y100 - 7;
1206  return 0.20476320420324610618e-1 + (0.30352843012898665856e-2 + (0.16617609387003727409e-4 + (0.13525429711163116103e-6 + (0.12515095552507169013e-8 + (0.13235687543603382345e-10 + 0.15326595042666666667e-12 * t) * t) * t) * t) * t) * t;
1207  }
1208  case 4 : {
1209  double t = 2*y100 - 9;
1210  return 0.26614461952489004566e-1 + (0.31034189276234947088e-2 + (0.17460268109986214274e-4 + (0.14582130824485709573e-6 + (0.13935959083809746345e-8 + (0.15249438072998932900e-10 + 0.18344741882133333333e-12 * t) * t) * t) * t) * t) * t;
1211  }
1212  case 5 : {
1213  double t = 2*y100 - 11;
1214  return 0.32892330248093586215e-1 + (0.31750557067975068584e-2 + (0.18369907582308672632e-4 + (0.15761063702089457882e-6 + (0.15577638230480894382e-8 + (0.17663868462699097951e-10 + (0.22126732680711111111e-12 + 0.30273474177737853668e-14 * t) * t) * t) * t) * t) * t) * t;
1215  }
1216  case 6 : {
1217  double t = 2*y100 - 13;
1218  return 0.39317207681134336024e-1 + (0.32504779701937539333e-2 + (0.19354426046513400534e-4 + (0.17081646971321290539e-6 + (0.17485733959327106250e-8 + (0.20593687304921961410e-10 + (0.26917401949155555556e-12 + 0.38562123837725712270e-14 * t) * t) * t) * t) * t) * t) * t;
1219  }
1220  case 7 : {
1221  double t = 2*y100 - 15;
1222  return 0.45896976511367738235e-1 + (0.33300031273110976165e-2 + (0.20423005398039037313e-4 + (0.18567412470376467303e-6 + (0.19718038363586588213e-8 + (0.24175006536781219807e-10 + (0.33059982791466666666e-12 + 0.49756574284439426165e-14 * t) * t) * t) * t) * t) * t) * t;
1223  }
1224  case 8 : {
1225  double t = 2*y100 - 17;
1226  return 0.52640192524848962855e-1 + (0.34139883358846720806e-2 + (0.21586390240603337337e-4 + (0.20247136501568904646e-6 + (0.22348696948197102935e-8 + (0.28597516301950162548e-10 + (0.41045502119111111110e-12 + 0.65151614515238361946e-14 * t) * t) * t) * t) * t) * t) * t;
1227  }
1228  case 9 : {
1229  double t = 2*y100 - 19;
1230  return 0.59556171228656770456e-1 + (0.35028374386648914444e-2 + (0.22857246150998562824e-4 + (0.22156372146525190679e-6 + (0.25474171590893813583e-8 + (0.34122390890697400584e-10 + (0.51593189879111111110e-12 + 0.86775076853908006938e-14 * t) * t) * t) * t) * t) * t) * t;
1231  }
1232  case 10 : {
1233  double t = 2*y100 - 21;
1234  return 0.66655089485108212551e-1 + (0.35970095381271285568e-2 + (0.24250626164318672928e-4 + (0.24339561521785040536e-6 + (0.29221990406518411415e-8 + (0.41117013527967776467e-10 + (0.65786450716444444445e-12 + 0.11791885745450623331e-13 * t) * t) * t) * t) * t) * t) * t;
1235  }
1236  case 11 : {
1237  double t = 2*y100 - 23;
1238  return 0.73948106345519174661e-1 + (0.36970297216569341748e-2 + (0.25784588137312868792e-4 + (0.26853012002366752770e-6 + (0.33763958861206729592e-8 + (0.50111549981376976397e-10 + (0.85313857496888888890e-12 + 0.16417079927706899860e-13 * t) * t) * t) * t) * t) * t) * t;
1239  }
1240  case 12 : {
1241  double t = 2*y100 - 25;
1242  return 0.81447508065002963203e-1 + (0.38035026606492705117e-2 + (0.27481027572231851896e-4 + (0.29769200731832331364e-6 + (0.39336816287457655076e-8 + (0.61895471132038157624e-10 + (0.11292303213511111111e-11 + 0.23558532213703884304e-13 * t) * t) * t) * t) * t) * t) * t;
1243  }
1244  case 13 : {
1245  double t = 2*y100 - 27;
1246  return 0.89166884027582716628e-1 + (0.39171301322438946014e-2 + (0.29366827260422311668e-4 + (0.33183204390350724895e-6 + (0.46276006281647330524e-8 + (0.77692631378169813324e-10 + (0.15335153258844444444e-11 + 0.35183103415916026911e-13 * t) * t) * t) * t) * t) * t) * t;
1247  }
1248  case 14 : {
1249  double t = 2*y100 - 29;
1250  return 0.97121342888032322019e-1 + (0.40387340353207909514e-2 + (0.31475490395950776930e-4 + (0.37222714227125135042e-6 + (0.55074373178613809996e-8 + (0.99509175283990337944e-10 + (0.21552645758222222222e-11 + 0.55728651431872687605e-13 * t) * t) * t) * t) * t) * t) * t;
1251  }
1252  case 15 : {
1253  double t = 2*y100 - 31;
1254  return 0.10532778218603311137e0 + (0.41692873614065380607e-2 + (0.33849549774889456984e-4 + (0.42064596193692630143e-6 + (0.66494579697622432987e-8 + (0.13094103581931802337e-9 + (0.31896187409777777778e-11 + 0.97271974184476560742e-13 * t) * t) * t) * t) * t) * t) * t;
1255  }
1256  case 16 : {
1257  double t = 2*y100 - 33;
1258  return 0.11380523107427108222e0 + (0.43099572287871821013e-2 + (0.36544324341565929930e-4 + (0.47965044028581857764e-6 + (0.81819034238463698796e-8 + (0.17934133239549647357e-9 + (0.50956666166186293627e-11 + (0.18850487318190638010e-12 + 0.79697813173519853340e-14 * t) * t) * t) * t) * t) * t) * t) * t;
1259  }
1260  case 17 : {
1261  double t = 2*y100 - 35;
1262  return 0.12257529703447467345e0 + (0.44621675710026986366e-2 + (0.39634304721292440285e-4 + (0.55321553769873381819e-6 + (0.10343619428848520870e-7 + (0.26033830170470368088e-9 + (0.87743837749108025357e-11 + (0.34427092430230063401e-12 + 0.10205506615709843189e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1263  }
1264  case 18 : {
1265  double t = 2*y100 - 37;
1266  return 0.13166276955656699478e0 + (0.46276970481783001803e-2 + (0.43225026380496399310e-4 + (0.64799164020016902656e-6 + (0.13580082794704641782e-7 + (0.39839800853954313927e-9 + (0.14431142411840000000e-10 + 0.42193457308830027541e-12 * t) * t) * t) * t) * t) * t) * t;
1267  }
1268  case 19 : {
1269  double t = 2*y100 - 39;
1270  return 0.14109647869803356475e0 + (0.48088424418545347758e-2 + (0.47474504753352150205e-4 + (0.77509866468724360352e-6 + (0.18536851570794291724e-7 + (0.60146623257887570439e-9 + (0.18533978397305276318e-10 + (0.41033845938901048380e-13 - 0.46160680279304825485e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1271  }
1272  case 20 : {
1273  double t = 2*y100 - 41;
1274  return 0.15091057940548936603e0 + (0.50086864672004685703e-2 + (0.52622482832192230762e-4 + (0.95034664722040355212e-6 + (0.25614261331144718769e-7 + (0.80183196716888606252e-9 + (0.12282524750534352272e-10 + (-0.10531774117332273617e-11 - 0.86157181395039646412e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1275  }
1276  case 21 : {
1277  double t = 2*y100 - 43;
1278  return 0.16114648116017010770e0 + (0.52314661581655369795e-2 + (0.59005534545908331315e-4 + (0.11885518333915387760e-5 + (0.33975801443239949256e-7 + (0.82111547144080388610e-9 + (-0.12357674017312854138e-10 + (-0.24355112256914479176e-11 - 0.75155506863572930844e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1279  }
1280  case 22 : {
1281  double t = 2*y100 - 45;
1282  return 0.17185551279680451144e0 + (0.54829002967599420860e-2 + (0.67013226658738082118e-4 + (0.14897400671425088807e-5 + (0.40690283917126153701e-7 + (0.44060872913473778318e-9 + (-0.52641873433280000000e-10 - 0.30940587864543343124e-11 * t) * t) * t) * t) * t) * t) * t;
1283  }
1284  case 23 : {
1285  double t = 2*y100 - 47;
1286  return 0.18310194559815257381e0 + (0.57701559375966953174e-2 + (0.76948789401735193483e-4 + (0.18227569842290822512e-5 + (0.41092208344387212276e-7 + (-0.44009499965694442143e-9 + (-0.92195414685628803451e-10 + (-0.22657389705721753299e-11 + 0.10004784908106839254e-12 * t) * t) * t) * t) * t) * t) * t) * t;
1287  }
1288  case 24 : {
1289  double t = 2*y100 - 49;
1290  return 0.19496527191546630345e0 + (0.61010853144364724856e-2 + (0.88812881056342004864e-4 + (0.21180686746360261031e-5 + (0.30652145555130049203e-7 + (-0.16841328574105890409e-8 + (-0.11008129460612823934e-9 + (-0.12180794204544515779e-12 + 0.15703325634590334097e-12 * t) * t) * t) * t) * t) * t) * t) * t;
1291  }
1292  case 25 : {
1293  double t = 2*y100 - 51;
1294  return 0.20754006813966575720e0 + (0.64825787724922073908e-2 + (0.10209599627522311893e-3 + (0.22785233392557600468e-5 + (0.73495224449907568402e-8 + (-0.29442705974150112783e-8 + (-0.94082603434315016546e-10 + (0.23609990400179321267e-11 + 0.14141908654269023788e-12 * t) * t) * t) * t) * t) * t) * t) * t;
1295  }
1296  case 26 : {
1297  double t = 2*y100 - 53;
1298  return 0.22093185554845172146e0 + (0.69182878150187964499e-2 + (0.11568723331156335712e-3 + (0.22060577946323627739e-5 + (-0.26929730679360840096e-7 + (-0.38176506152362058013e-8 + (-0.47399503861054459243e-10 + (0.40953700187172127264e-11 + 0.69157730376118511127e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1299  }
1300  case 27 : {
1301  double t = 2*y100 - 55;
1302  return 0.23524827304057813918e0 + (0.74063350762008734520e-2 + (0.12796333874615790348e-3 + (0.18327267316171054273e-5 + (-0.66742910737957100098e-7 + (-0.40204740975496797870e-8 + (0.14515984139495745330e-10 + (0.44921608954536047975e-11 - 0.18583341338983776219e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1303  }
1304  case 28 : {
1305  double t = 2*y100 - 57;
1306  return 0.25058626331812744775e0 + (0.79377285151602061328e-2 + (0.13704268650417478346e-3 + (0.11427511739544695861e-5 + (-0.10485442447768377485e-6 + (-0.34850364756499369763e-8 + (0.72656453829502179208e-10 + (0.36195460197779299406e-11 - 0.84882136022200714710e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1307  }
1308  case 29 : {
1309  double t = 2*y100 - 59;
1310  return 0.26701724900280689785e0 + (0.84959936119625864274e-2 + (0.14112359443938883232e-3 + (0.17800427288596909634e-6 + (-0.13443492107643109071e-6 + (-0.23512456315677680293e-8 + (0.11245846264695936769e-9 + (0.19850501334649565404e-11 - 0.11284666134635050832e-12 * t) * t) * t) * t) * t) * t) * t) * t;
1311  }
1312  case 30 : {
1313  double t = 2*y100 - 61;
1314  return 0.28457293586253654144e0 + (0.90581563892650431899e-2 + (0.13880520331140646738e-3 + (-0.97262302362522896157e-6 + (-0.15077100040254187366e-6 + (-0.88574317464577116689e-9 + (0.12760311125637474581e-9 + (0.20155151018282695055e-12 - 0.10514169375181734921e-12 * t) * t) * t) * t) * t) * t) * t) * t;
1315  }
1316  case 31 : {
1317  double t = 2*y100 - 63;
1318  return 0.30323425595617385705e0 + (0.95968346790597422934e-2 + (0.12931067776725883939e-3 + (-0.21938741702795543986e-5 + (-0.15202888584907373963e-6 + (0.61788350541116331411e-9 + (0.11957835742791248256e-9 + (-0.12598179834007710908e-11 - 0.75151817129574614194e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1319  }
1320  case 32 : {
1321  double t = 2*y100 - 65;
1322  return 0.32292521181517384379e0 + (0.10082957727001199408e-1 + (0.11257589426154962226e-3 + (-0.33670890319327881129e-5 + (-0.13910529040004008158e-6 + (0.19170714373047512945e-8 + (0.94840222377720494290e-10 + (-0.21650018351795353201e-11 - 0.37875211678024922689e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1323  }
1324  case 33 : {
1325  double t = 2*y100 - 67;
1326  return 0.34351233557911753862e0 + (0.10488575435572745309e-1 + (0.89209444197248726614e-4 + (-0.43893459576483345364e-5 + (-0.11488595830450424419e-6 + (0.28599494117122464806e-8 + (0.61537542799857777779e-10 - 0.24935749227658002212e-11 * t) * t) * t) * t) * t) * t) * t;
1327  }
1328  case 34 : {
1329  double t = 2*y100 - 69;
1330  return 0.36480946642143669093e0 + (0.10789304203431861366e-1 + (0.60357993745283076834e-4 + (-0.51855862174130669389e-5 + (-0.83291664087289801313e-7 + (0.33898011178582671546e-8 + (0.27082948188277716482e-10 + (-0.23603379397408694974e-11 + 0.19328087692252869842e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1331  }
1332  case 35 : {
1333  double t = 2*y100 - 71;
1334  return 0.38658679935694939199e0 + (0.10966119158288804999e-1 + (0.27521612041849561426e-4 + (-0.57132774537670953638e-5 + (-0.48404772799207914899e-7 + (0.35268354132474570493e-8 + (-0.32383477652514618094e-11 + (-0.19334202915190442501e-11 + 0.32333189861286460270e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1335  }
1336  case 36 : {
1337  double t = 2*y100 - 73;
1338  return 0.40858275583808707870e0 + (0.11006378016848466550e-1 + (-0.76396376685213286033e-5 + (-0.59609835484245791439e-5 + (-0.13834610033859313213e-7 + (0.33406952974861448790e-8 + (-0.26474915974296612559e-10 + (-0.13750229270354351983e-11 + 0.36169366979417390637e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1339  }
1340  case 37 : {
1341  double t = 2*y100 - 75;
1342  return 0.43051714914006682977e0 + (0.10904106549500816155e-1 + (-0.43477527256787216909e-4 + (-0.59429739547798343948e-5 + (0.17639200194091885949e-7 + (0.29235991689639918688e-8 + (-0.41718791216277812879e-10 + (-0.81023337739508049606e-12 + 0.33618915934461994428e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1343  }
1344  case 38 : {
1345  double t = 2*y100 - 77;
1346  return 0.45210428135559607406e0 + (0.10659670756384400554e-1 + (-0.78488639913256978087e-4 + (-0.56919860886214735936e-5 + (0.44181850467477733407e-7 + (0.23694306174312688151e-8 + (-0.49492621596685443247e-10 + (-0.31827275712126287222e-12 + 0.27494438742721623654e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1347  }
1348  case 39 : {
1349  double t = 2*y100 - 79;
1350  return 0.47306491195005224077e0 + (0.10279006119745977570e-1 + (-0.11140268171830478306e-3 + (-0.52518035247451432069e-5 + (0.64846898158889479518e-7 + (0.17603624837787337662e-8 + (-0.51129481592926104316e-10 + (0.62674584974141049511e-13 + 0.20055478560829935356e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1351  }
1352  case 40 : {
1353  double t = 2*y100 - 81;
1354  return 0.49313638965719857647e0 + (0.97725799114772017662e-2 + (-0.14122854267291533334e-3 + (-0.46707252568834951907e-5 + (0.79421347979319449524e-7 + (0.11603027184324708643e-8 + (-0.48269605844397175946e-10 + (0.32477251431748571219e-12 + 0.12831052634143527985e-13 * t) * t) * t) * t) * t) * t) * t) * t;
1355  }
1356  case 41 : {
1357  double t = 2*y100 - 83;
1358  return 0.51208057433416004042e0 + (0.91542422354009224951e-2 + (-0.16726530230228647275e-3 + (-0.39964621752527649409e-5 + (0.88232252903213171454e-7 + (0.61343113364949928501e-9 + (-0.42516755603130443051e-10 + (0.47910437172240209262e-12 + 0.66784341874437478953e-14 * t) * t) * t) * t) * t) * t) * t) * t;
1359  }
1360  case 42 : {
1361  double t = 2*y100 - 85;
1362  return 0.52968945458607484524e0 + (0.84400880445116786088e-2 + (-0.18908729783854258774e-3 + (-0.32725905467782951931e-5 + (0.91956190588652090659e-7 + (0.14593989152420122909e-9 + (-0.35239490687644444445e-10 + 0.54613829888448694898e-12 * t) * t) * t) * t) * t) * t) * t;
1363  }
1364  case 43 : {
1365  double t = 2*y100 - 87;
1366  return 0.54578857454330070965e0 + (0.76474155195880295311e-2 + (-0.20651230590808213884e-3 + (-0.25364339140543131706e-5 + (0.91455367999510681979e-7 + (-0.23061359005297528898e-9 + (-0.27512928625244444444e-10 + 0.54895806008493285579e-12 * t) * t) * t) * t) * t) * t) * t;
1367  }
1368  case 44 : {
1369  double t = 2*y100 - 89;
1370  return 0.56023851910298493910e0 + (0.67938321739997196804e-2 + (-0.21956066613331411760e-3 + (-0.18181127670443266395e-5 + (0.87650335075416845987e-7 + (-0.51548062050366615977e-9 + (-0.20068462174044444444e-10 + 0.50912654909758187264e-12 * t) * t) * t) * t) * t) * t) * t;
1371  }
1372  case 45 : {
1373  double t = 2*y100 - 91;
1374  return 0.57293478057455721150e0 + (0.58965321010394044087e-2 + (-0.22841145229276575597e-3 + (-0.11404605562013443659e-5 + (0.81430290992322326296e-7 + (-0.71512447242755357629e-9 + (-0.13372664928000000000e-10 + 0.44461498336689298148e-12 * t) * t) * t) * t) * t) * t) * t;
1375  }
1376  case 46 : {
1377  double t = 2*y100 - 93;
1378  return 0.58380635448407827360e0 + (0.49717469530842831182e-2 + (-0.23336001540009645365e-3 + (-0.51952064448608850822e-6 + (0.73596577815411080511e-7 + (-0.84020916763091566035e-9 + (-0.76700972702222222221e-11 + 0.36914462807972467044e-12 * t) * t) * t) * t) * t) * t) * t;
1379  }
1380  case 47 : {
1381  double t = 2*y100 - 95;
1382  return 0.59281340237769489597e0 + (0.40343592069379730568e-2 + (-0.23477963738658326185e-3 + (0.34615944987790224234e-7 + (0.64832803248395814574e-7 + (-0.90329163587627007971e-9 + (-0.30421940400000000000e-11 + 0.29237386653743536669e-12 * t) * t) * t) * t) * t) * t) * t;
1383  }
1384  case 48 : {
1385  double t = 2*y100 - 97;
1386  return 0.59994428743114271918e0 + (0.30976579788271744329e-2 + (-0.23308875765700082835e-3 + (0.51681681023846925160e-6 + (0.55694594264948268169e-7 + (-0.91719117313243464652e-9 + (0.53982743680000000000e-12 + 0.22050829296187771142e-12 * t) * t) * t) * t) * t) * t) * t;
1387  }
1388  case 49 : {
1389  double t = 2*y100 - 99;
1390  return 0.60521224471819875444e0 + (0.21732138012345456060e-2 + (-0.22872428969625997456e-3 + (0.92588959922653404233e-6 + (0.46612665806531930684e-7 + (-0.89393722514414153351e-9 + (0.31718550353777777778e-11 + 0.15705458816080549117e-12 * t) * t) * t) * t) * t) * t) * t;
1391  }
1392  case 50 : {
1393  double t = 2*y100 - 101;
1394  return 0.60865189969791123620e0 + (0.12708480848877451719e-2 + (-0.22212090111534847166e-3 + (0.12636236031532793467e-5 + (0.37904037100232937574e-7 + (-0.84417089968101223519e-9 + (0.49843180828444444445e-11 + 0.10355439441049048273e-12 * t) * t) * t) * t) * t) * t) * t;
1395  }
1396  case 51 : {
1397  double t = 2*y100 - 103;
1398  return 0.61031580103499200191e0 + (0.39867436055861038223e-3 + (-0.21369573439579869291e-3 + (0.15339402129026183670e-5 + (0.29787479206646594442e-7 + (-0.77687792914228632974e-9 + (0.61192452741333333334e-11 + 0.60216691829459295780e-13 * t) * t) * t) * t) * t) * t) * t;
1399  }
1400  case 52 : {
1401  double t = 2*y100 - 105;
1402  return 0.61027109047879835868e0 + (-0.43680904508059878254e-3 + (-0.20383783788303894442e-3 + (0.17421743090883439959e-5 + (0.22400425572175715576e-7 + (-0.69934719320045128997e-9 + (0.67152759655111111110e-11 + 0.26419960042578359995e-13 * t) * t) * t) * t) * t) * t) * t;
1403  }
1404  case 53 : {
1405  double t = 2*y100 - 107;
1406  return 0.60859639489217430521e0 + (-0.12305921390962936873e-2 + (-0.19290150253894682629e-3 + (0.18944904654478310128e-5 + (0.15815530398618149110e-7 + (-0.61726850580964876070e-9 + 0.68987888999111111110e-11 * t) * t) * t) * t) * t) * t;
1407  }
1408  case 54 : {
1409  double t = 2*y100 - 109;
1410  return 0.60537899426486075181e0 + (-0.19790062241395705751e-2 + (-0.18120271393047062253e-3 + (0.19974264162313241405e-5 + (0.10055795094298172492e-7 + (-0.53491997919318263593e-9 + (0.67794550295111111110e-11 - 0.17059208095741511603e-13 * t) * t) * t) * t) * t) * t) * t;
1411  }
1412  case 55 : {
1413  double t = 2*y100 - 111;
1414  return 0.60071229457904110537e0 + (-0.26795676776166354354e-2 + (-0.16901799553627508781e-3 + (0.20575498324332621581e-5 + (0.51077165074461745053e-8 + (-0.45536079828057221858e-9 + (0.64488005516444444445e-11 - 0.29311677573152766338e-13 * t) * t) * t) * t) * t) * t) * t;
1415  }
1416  case 56 : {
1417  double t = 2*y100 - 113;
1418  return 0.59469361520112714738e0 + (-0.33308208190600993470e-2 + (-0.15658501295912405679e-3 + (0.20812116912895417272e-5 + (0.93227468760614182021e-9 + (-0.38066673740116080415e-9 + (0.59806790359111111110e-11 - 0.36887077278950440597e-13 * t) * t) * t) * t) * t) * t) * t;
1419  }
1420  case 57 : {
1421  double t = 2*y100 - 115;
1422  return 0.58742228631775388268e0 + (-0.39321858196059227251e-2 + (-0.14410441141450122535e-3 + (0.20743790018404020716e-5 + (-0.25261903811221913762e-8 + (-0.31212416519526924318e-9 + (0.54328422462222222221e-11 - 0.40864152484979815972e-13 * t) * t) * t) * t) * t) * t) * t;
1423  }
1424  case 58 : {
1425  double t = 2*y100 - 117;
1426  return 0.57899804200033018447e0 + (-0.44838157005618913447e-2 + (-0.13174245966501437965e-3 + (0.20425306888294362674e-5 + (-0.53330296023875447782e-8 + (-0.25041289435539821014e-9 + (0.48490437205333333334e-11 - 0.42162206939169045177e-13 * t) * t) * t) * t) * t) * t) * t;
1427  }
1428  case 59 : {
1429  double t = 2*y100 - 119;
1430  return 0.56951968796931245974e0 + (-0.49864649488074868952e-2 + (-0.11963416583477567125e-3 + (0.19906021780991036425e-5 + (-0.75580140299436494248e-8 + (-0.19576060961919820491e-9 + (0.42613011928888888890e-11 - 0.41539443304115604377e-13 * t) * t) * t) * t) * t) * t) * t;
1431  }
1432  case 60 : {
1433  double t = 2*y100 - 121;
1434  return 0.55908401930063918964e0 + (-0.54413711036826877753e-2 + (-0.10788661102511914628e-3 + (0.19229663322982839331e-5 + (-0.92714731195118129616e-8 + (-0.14807038677197394186e-9 + (0.36920870298666666666e-11 - 0.39603726688419162617e-13 * t) * t) * t) * t) * t) * t) * t;
1435  }
1436  case 61 : {
1437  double t = 2*y100 - 123;
1438  return 0.54778496152925675315e0 + (-0.58501497933213396670e-2 + (-0.96582314317855227421e-4 + (0.18434405235069270228e-5 + (-0.10541580254317078711e-7 + (-0.10702303407788943498e-9 + (0.31563175582222222222e-11 - 0.36829748079110481422e-13 * t) * t) * t) * t) * t) * t) * t;
1439  }
1440  case 62 : {
1441  double t = 2*y100 - 125;
1442  return 0.53571290831682823999e0 + (-0.62147030670760791791e-2 + (-0.85782497917111760790e-4 + (0.17553116363443470478e-5 + (-0.11432547349815541084e-7 + (-0.72157091369041330520e-10 + (0.26630811607111111111e-11 - 0.33578660425893164084e-13 * t) * t) * t) * t) * t) * t) * t;
1443  }
1444  case 63 : {
1445  double t = 2*y100 - 127;
1446  return 0.52295422962048434978e0 + (-0.65371404367776320720e-2 + (-0.75530164941473343780e-4 + (0.16613725797181276790e-5 + (-0.12003521296598910761e-7 + (-0.42929753689181106171e-10 + (0.22170894940444444444e-11 - 0.30117697501065110505e-13 * t) * t) * t) * t) * t) * t) * t;
1447  }
1448  case 64 : {
1449  double t = 2*y100 - 129;
1450  return 0.50959092577577886140e0 + (-0.68197117603118591766e-2 + (-0.65852936198953623307e-4 + (0.15639654113906716939e-5 + (-0.12308007991056524902e-7 + (-0.18761997536910939570e-10 + (0.18198628922666666667e-11 - 0.26638355362285200932e-13 * t) * t) * t) * t) * t) * t) * t;
1451  }
1452  case 65 : {
1453  double t = 2*y100 - 131;
1454  return 0.49570040481823167970e0 + (-0.70647509397614398066e-2 + (-0.56765617728962588218e-4 + (0.14650274449141448497e-5 + (-0.12393681471984051132e-7 + (0.92904351801168955424e-12 + (0.14706755960177777778e-11 - 0.23272455351266325318e-13 * t) * t) * t) * t) * t) * t) * t;
1455  }
1456  case 66 : {
1457  double t = 2*y100 - 133;
1458  return 0.48135536250935238066e0 + (-0.72746293327402359783e-2 + (-0.48272489495730030780e-4 + (0.13661377309113939689e-5 + (-0.12302464447599382189e-7 + (0.16707760028737074907e-10 + (0.11672928324444444444e-11 - 0.20105801424709924499e-13 * t) * t) * t) * t) * t) * t) * t;
1459  }
1460  case 67 : {
1461  double t = 2*y100 - 135;
1462  return 0.46662374675511439448e0 + (-0.74517177649528487002e-2 + (-0.40369318744279128718e-4 + (0.12685621118898535407e-5 + (-0.12070791463315156250e-7 + (0.29105507892605823871e-10 + (0.90653314645333333334e-12 - 0.17189503312102982646e-13 * t) * t) * t) * t) * t) * t) * t;
1463  }
1464  case 68 : {
1465  double t = 2*y100 - 137;
1466  return 0.45156879030168268778e0 + (-0.75983560650033817497e-2 + (-0.33045110380705139759e-4 + (0.11732956732035040896e-5 + (-0.11729986947158201869e-7 + (0.38611905704166441308e-10 + (0.68468768305777777779e-12 - 0.14549134330396754575e-13 * t) * t) * t) * t) * t) * t) * t;
1467  }
1468  case 69 : {
1469  double t = 2*y100 - 139;
1470  return 0.43624909769330896904e0 + (-0.77168291040309554679e-2 + (-0.26283612321339907756e-4 + (0.10811018836893550820e-5 + (-0.11306707563739851552e-7 + (0.45670446788529607380e-10 + (0.49782492549333333334e-12 - 0.12191983967561779442e-13 * t) * t) * t) * t) * t) * t) * t;
1471  }
1472  case 70 : {
1473  double t = 2*y100 - 141;
1474  return 0.42071877443548481181e0 + (-0.78093484015052730097e-2 + (-0.20064596897224934705e-4 + (0.99254806680671890766e-6 + (-0.10823412088884741451e-7 + (0.50677203326904716247e-10 + (0.34200547594666666666e-12 - 0.10112698698356194618e-13 * t) * t) * t) * t) * t) * t) * t;
1475  }
1476  case 71 : {
1477  double t = 2*y100 - 143;
1478  return 0.40502758809710844280e0 + (-0.78780384460872937555e-2 + (-0.14364940764532853112e-4 + (0.90803709228265217384e-6 + (-0.10298832847014466907e-7 + (0.53981671221969478551e-10 + (0.21342751381333333333e-12 - 0.82975901848387729274e-14 * t) * t) * t) * t) * t) * t) * t;
1479  }
1480  case 72 : {
1481  double t = 2*y100 - 145;
1482  return 0.38922115269731446690e0 + (-0.79249269708242064120e-2 + (-0.91595258799106970453e-5 + (0.82783535102217576495e-6 + (-0.97484311059617744437e-8 + (0.55889029041660225629e-10 + (0.10851981336888888889e-12 - 0.67278553237853459757e-14 * t) * t) * t) * t) * t) * t) * t;
1483  }
1484  case 73 : {
1485  double t = 2*y100 - 147;
1486  return 0.37334112915460307335e0 + (-0.79519385109223148791e-2 + (-0.44219833548840469752e-5 + (0.75209719038240314732e-6 + (-0.91848251458553190451e-8 + (0.56663266668051433844e-10 + (0.23995894257777777778e-13 - 0.53819475285389344313e-14 * t) * t) * t) * t) * t) * t) * t;
1487  }
1488  case 74 : {
1489  double t = 2*y100 - 149;
1490  return 0.35742543583374223085e0 + (-0.79608906571527956177e-2 + (-0.12530071050975781198e-6 + (0.68088605744900552505e-6 + (-0.86181844090844164075e-8 + (0.56530784203816176153e-10 + (-0.43120012248888888890e-13 - 0.42372603392496813810e-14 * t) * t) * t) * t) * t) * t) * t;
1491  }
1492  case 75 : {
1493  double t = 2*y100 - 151;
1494  return 0.34150846431979618536e0 + (-0.79534924968773806029e-2 + (0.37576885610891515813e-5 + (0.61419263633090524326e-6 + (-0.80565865409945960125e-8 + (0.55684175248749269411e-10 + (-0.95486860764444444445e-13 - 0.32712946432984510595e-14 * t) * t) * t) * t) * t) * t) * t;
1495  }
1496  case 76 : {
1497  double t = 2*y100 - 153;
1498  return 0.32562129649136346824e0 + (-0.79313448067948884309e-2 + (0.72539159933545300034e-5 + (0.55195028297415503083e-6 + (-0.75063365335570475258e-8 + (0.54281686749699595941e-10 - 0.13545424295111111111e-12 * t) * t) * t) * t) * t) * t;
1499  }
1500  case 77 : {
1501  double t = 2*y100 - 155;
1502  return 0.30979191977078391864e0 + (-0.78959416264207333695e-2 + (0.10389774377677210794e-4 + (0.49404804463196316464e-6 + (-0.69722488229411164685e-8 + (0.52469254655951393842e-10 - 0.16507860650666666667e-12 * t) * t) * t) * t) * t) * t;
1503  }
1504  case 78 : {
1505  double t = 2*y100 - 157;
1506  return 0.29404543811214459904e0 + (-0.78486728990364155356e-2 + (0.13190885683106990459e-4 + (0.44034158861387909694e-6 + (-0.64578942561562616481e-8 + (0.50354306498006928984e-10 - 0.18614473550222222222e-12 * t) * t) * t) * t) * t) * t;
1507  }
1508  case 79 : {
1509  double t = 2*y100 - 159;
1510  return 0.27840427686253660515e0 + (-0.77908279176252742013e-2 + (0.15681928798708548349e-4 + (0.39066226205099807573e-6 + (-0.59658144820660420814e-8 + (0.48030086420373141763e-10 - 0.20018995173333333333e-12 * t) * t) * t) * t) * t) * t;
1511  }
1512  case 80 : {
1513  double t = 2*y100 - 161;
1514  return 0.26288838011163800908e0 + (-0.77235993576119469018e-2 + (0.17886516796198660969e-4 + (0.34482457073472497720e-6 + (-0.54977066551955420066e-8 + (0.45572749379147269213e-10 - 0.20852924954666666667e-12 * t) * t) * t) * t) * t) * t;
1515  }
1516  case 81 : {
1517  double t = 2*y100 - 163;
1518  return 0.24751539954181029717e0 + (-0.76480877165290370975e-2 + (0.19827114835033977049e-4 + (0.30263228619976332110e-6 + (-0.50545814570120129947e-8 + (0.43043879374212005966e-10 - 0.21228012028444444444e-12 * t) * t) * t) * t) * t) * t;
1519  }
1520  case 82 : {
1521  double t = 2*y100 - 165;
1522  return 0.23230087411688914593e0 + (-0.75653060136384041587e-2 + (0.21524991113020016415e-4 + (0.26388338542539382413e-6 + (-0.46368974069671446622e-8 + (0.40492715758206515307e-10 - 0.21238627815111111111e-12 * t) * t) * t) * t) * t) * t;
1523  }
1524  case 83 : {
1525  double t = 2*y100 - 167;
1526  return 0.21725840021297341931e0 + (-0.74761846305979730439e-2 + (0.23000194404129495243e-4 + (0.22837400135642906796e-6 + (-0.42446743058417541277e-8 + (0.37958104071765923728e-10 - 0.20963978568888888889e-12 * t) * t) * t) * t) * t) * t;
1527  }
1528  case 84 : {
1529  double t = 2*y100 - 169;
1530  return 0.20239979200788191491e0 + (-0.73815761980493466516e-2 + (0.24271552727631854013e-4 + (0.19590154043390012843e-6 + (-0.38775884642456551753e-8 + (0.35470192372162901168e-10 - 0.20470131678222222222e-12 * t) * t) * t) * t) * t) * t;
1531  }
1532  case 85 : {
1533  double t = 2*y100 - 171;
1534  return 0.18773523211558098962e0 + (-0.72822604530339834448e-2 + (0.25356688567841293697e-4 + (0.16626710297744290016e-6 + (-0.35350521468015310830e-8 + (0.33051896213898864306e-10 - 0.19811844544000000000e-12 * t) * t) * t) * t) * t) * t;
1535  }
1536  case 86 : {
1537  double t = 2*y100 - 173;
1538  return 0.17327341258479649442e0 + (-0.71789490089142761950e-2 + (0.26272046822383820476e-4 + (0.13927732375657362345e-6 + (-0.32162794266956859603e-8 + (0.30720156036105652035e-10 - 0.19034196304000000000e-12 * t) * t) * t) * t) * t) * t;
1539  }
1540  case 87 : {
1541  double t = 2*y100 - 175;
1542  return 0.15902166648328672043e0 + (-0.70722899934245504034e-2 + (0.27032932310132226025e-4 + (0.11474573347816568279e-6 + (-0.29203404091754665063e-8 + (0.28487010262547971859e-10 - 0.18174029063111111111e-12 * t) * t) * t) * t) * t) * t;
1543  }
1544  case 88 : {
1545  double t = 2*y100 - 177;
1546  return 0.14498609036610283865e0 + (-0.69628725220045029273e-2 + (0.27653554229160596221e-4 + (0.92493727167393036470e-7 + (-0.26462055548683583849e-8 + (0.26360506250989943739e-10 - 0.17261211260444444444e-12 * t) * t) * t) * t) * t) * t;
1547  }
1548  case 89 : {
1549  double t = 2*y100 - 179;
1550  return 0.13117165798208050667e0 + (-0.68512309830281084723e-2 + (0.28147075431133863774e-4 + (0.72351212437979583441e-7 + (-0.23927816200314358570e-8 + (0.24345469651209833155e-10 - 0.16319736960000000000e-12 * t) * t) * t) * t) * t) * t;
1551  }
1552  case 90 : {
1553  double t = 2*y100 - 181;
1554  return 0.11758232561160626306e0 + (-0.67378491192463392927e-2 + (0.28525664781722907847e-4 + (0.54156999310046790024e-7 + (-0.21589405340123827823e-8 + (0.22444150951727334619e-10 - 0.15368675584000000000e-12 * t) * t) * t) * t) * t) * t;
1555  }
1556  case 91 : {
1557  double t = 2*y100 - 183;
1558  return 0.10422112945361673560e0 + (-0.66231638959845581564e-2 + (0.28800551216363918088e-4 + (0.37758983397952149613e-7 + (-0.19435423557038933431e-8 + (0.20656766125421362458e-10 - 0.14422990012444444444e-12 * t) * t) * t) * t) * t) * t;
1559  }
1560  case 92 : {
1561  double t = 2*y100 - 185;
1562  return 0.91090275493541084785e-1 + (-0.65075691516115160062e-2 + (0.28982078385527224867e-4 + (0.23014165807643012781e-7 + (-0.17454532910249875958e-8 + (0.18981946442680092373e-10 - 0.13494234691555555556e-12 * t) * t) * t) * t) * t) * t;
1563  }
1564  case 93 : {
1565  double t = 2*y100 - 187;
1566  return 0.78191222288771379358e-1 + (-0.63914190297303976434e-2 + (0.29079759021299682675e-4 + (0.97885458059415717014e-8 + (-0.15635596116134296819e-8 + (0.17417110744051331974e-10 - 0.12591151763555555556e-12 * t) * t) * t) * t) * t) * t;
1567  }
1568  case 94 : {
1569  double t = 2*y100 - 189;
1570  return 0.65524757106147402224e-1 + (-0.62750311956082444159e-2 + (0.29102328354323449795e-4 + (-0.20430838882727954582e-8 + (-0.13967781903855367270e-8 + (0.15958771833747057569e-10 - 0.11720175765333333333e-12 * t) * t) * t) * t) * t) * t;
1571  }
1572  case 95 : {
1573  double t = 2*y100 - 191;
1574  return 0.53091065838453612773e-1 + (-0.61586898417077043662e-2 + (0.29057796072960100710e-4 + (-0.12597414620517987536e-7 + (-0.12440642607426861943e-8 + (0.14602787128447932137e-10 - 0.10885859114666666667e-12 * t) * t) * t) * t) * t) * t;
1575  }
1576  case 96 : {
1577  double t = 2*y100 - 193;
1578  return 0.40889797115352738582e-1 + (-0.60426484889413678200e-2 + (0.28953496450191694606e-4 + (-0.21982952021823718400e-7 + (-0.11044169117553026211e-8 + (0.13344562332430552171e-10 - 0.10091231402844444444e-12 * t) * t) * t) * t) * t) * t;
1579  }
1580  case 97: case 98:
1581  case 99: case 100 : { // use Taylor expansion for small x (|x| <= 0.0309...)
1582  // (2/sqrt(pi)) * (x - 2/3 x^3 + 4/15 x^5 - 8/105 x^7 + 16/945 x^9)
1583  double x2 = x*x;
1584  return x * (1.1283791670955125739
1585  - x2 * (0.75225277806367504925
1586  - x2 * (0.30090111122547001970
1587  - x2 * (0.085971746064420005629
1588  - x2 * 0.016931216931216931217))));
1589  }
1590  }
1591  // Since 0 <= y100 < 101, this is only reached if x is NaN,
1592  // in which case we should return NaN
1593  return NaN;
1594 }
1595 
1596 double w_im(double x) {
1597  if ( x >= 0 ) {
1598  if ( x > 45 ) { // continued-fraction expansion is faster
1599  const double ispi = 0.56418958354775628694807945156; // 1 / sqrt(pi)
1600  if ( x > 5e7 ) { // 1-term expansion, important to avoid overflow
1601  return ispi / x;
1602  }
1603  // 5-term expansion (rely on compiler for CSE), simplified from: ispi / (x-0.5/(x-1/(x-1.5/(x-2/x))))
1604  return ispi*((x*x) * (x*x-4.5) + 2) / (x * ((x*x) * (x*x-5) + 3.75));
1605  }
1606  return w_im_y100(100/(1+x), x);
1607  } else { // = -w_im(-x)
1608  if ( x < -45 ) { // continued-fraction expansion is faster
1609  const double ispi = 0.56418958354775628694807945156; // 1 / sqrt(pi)
1610  if ( x < -5e7 ) { // 1-term expansion, important to avoid overflow
1611  return ispi / x;
1612  }
1613  /* 5-term expansion (rely on compiler for CSE), simplified from:
1614  ispi / (x-0.5/(x-1/(x-1.5/(x-2/x)))) */
1615  return ispi*((x*x) * (x*x-4.5) + 2) / (x * ((x*x) * (x*x-5) + 3.75));
1616  }
1617  return -w_im_y100(100/(1-x), -x);
1618  }
1619 }
1620 
1621 
1622 } // namespace statistics
1623 } // namespace numeric
static double sinh_taylor(double x)
Definition: functions.cc:440
double Dawson(double x)
Definition: functions.cc:297
static double sinc(double x, double sinx)
Definition: functions.cc:435
platform::Size Size
Definition: types.hh:42
double errf(double x)
Definition: functions.cc:131
cmplx w(cmplx z, double relerr)
Definition: functions.cc:470
#define creal(z)
Definition: functions.cc:29
def x
#define cimag(z)
Definition: functions.cc:30
cmplx errfcx(cmplx z, double relerr)
Definition: functions.cc:126
double log(double x, double base)
Computes log(x) in the given base.
Definition: util.hh:41
bool is_inf(platform::Real const &val)
Definition: numbers.hh:66
static double sqr(double x)
Definition: functions.cc:445
Brief utility classes for numeric usage.
#define C(a, b)
Definition: functions.cc:27
numeric::Real kl_divergence(utility::vector1< numeric::Real > const &prior, utility::vector1< numeric::Real > const &posterior)
Returns the Kullback-Leibler divergence (aka relative entropy) between two discrete probability distr...
Definition: functions.cc:38
cmplx errfi(cmplx z, double relerr)
Definition: functions.cc:240
tuple p
Definition: docking.py:9
static double erfcx_y100(double y100)
Definition: functions.cc:747
numeric::Real corrcoef_with_provided_mean_and_std_dev(utility::vector1< numeric::Real > const &vec1, numeric::Real m1, numeric::Real sd1, utility::vector1< numeric::Real > const &vec2, numeric::Real m2, numeric::Real sd2)
Definition: functions.cc:71
#define cexp(z)
Definition: functions.cc:28
def z
numeric::Real corrcoef(utility::vector1< numeric::Real > const &vec1, utility::vector1< numeric::Real > const &vec2)
Definition: functions.cc:58
platform::Real copysign(platform::Real const &x, platform::Real const &y)
Definition: numbers.hh:76
std::complex< double > cmplx
Definition: functions.cc:26
std::vector with 1-based indexing
Definition: vector1.fwd.hh:44
numeric::Real cov(utility::vector1< numeric::Real > const &vec1, utility::vector1< numeric::Real > const &vec2)
Definition: functions.cc:91
bool is_nan(platform::Real const &val)
Definition: numbers.hh:61
static const unsigned int expa2n2_length
Definition: functions.cc:466
a collection of various functions to compute statistics. feel free to add your own ...
super::const_iterator const_iterator
Definition: vector1.hh:62
#define NaN
Definition: functions.cc:33
rosetta project type declarations. Should be kept updated with core/types.hh. This exists because num...
double w_im(double x)
Definition: functions.cc:1596
double Real
Definition: types.hh:39
T mean(Iterator first, Iterator last, T)
mean value of an input vector
Definition: functions.hh:34
static double w_im_y100(double y100, double x)
Definition: functions.cc:1190
T std_dev_with_provided_mean(Iterator first, Iterator last, T mean)
Definition: functions.hh:50
static const double expa2n2[]
Definition: functions.cc:449
vector1: std::vector with 1-based indexing
double errfc(double x)
Definition: functions.cc:252
#define pi
void statistics(std::string filename)
#define Inf
Definition: functions.cc:32
def y
numeric::Real cov_with_provided_mean(utility::vector1< numeric::Real > const &vec1, numeric::Real m1, utility::vector1< numeric::Real > const &vec2, numeric::Real m2)
Definition: functions.cc:101