Qucs-core  0.0.18
 All Data Structures Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
fspecial.cpp
Go to the documentation of this file.
1 /*
2  * fspecial.cpp - special functions implementation
3  *
4  * Copyright (C) 2006, 2008 Stefan Jahn <stefan@lkcc.org>
5  *
6  * This is free software; you can redistribute it and/or modify
7  * it under the terms of the GNU General Public License as published by
8  * the Free Software Foundation; either version 2, or (at your option)
9  * any later version.
10  *
11  * This software is distributed in the hope that it will be useful,
12  * but WITHOUT ANY WARRANTY; without even the implied warranty of
13  * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  * GNU General Public License for more details.
15  *
16  * You should have received a copy of the GNU General Public License
17  * along with this package; see the file COPYING. If not, write to
18  * the Free Software Foundation, Inc., 51 Franklin Street - Fifth Floor,
19  * Boston, MA 02110-1301, USA.
20  *
21  * $Id$
22  *
23  */
24 
25 #if HAVE_CONFIG_H
26 # include <config.h>
27 #endif
28 
29 #include <stdio.h>
30 #include <stdlib.h>
31 #include <string.h>
32 #include <cmath>
33 #include <float.h>
34 
35 #include "compat.h"
36 #include "constants.h"
37 #include "fspecial.h"
38 
39 /* The function computes the complete elliptic integral of first kind
40  K() and the second kind E() using the arithmetic-geometric mean
41  algorithm (AGM) by Abramowitz and Stegun. */
42 void fspecial::ellip_ke (nr_double_t arg, nr_double_t &k, nr_double_t &e) {
43  int iMax = 16;
44  if (arg == 1.0) {
45  k = NR_INF; // infinite
46  e = 0;
47  }
48  else if (std::isinf (arg) && arg < 0) {
49  k = 0;
50  e = NR_INF; // infinite
51  }
52  else {
53  nr_double_t a, b, c, f, s, fk = 1, fe = 1, t, da = arg;
54  int i;
55  if (arg < 0) {
56  fk = 1 / sqrt (1 - arg);
57  fe = sqrt (1 - arg);
58  da = -arg / (1 - arg);
59  }
60  a = 1;
61  b = sqrt (1 - da);
62  c = sqrt (da);
63  f = 0.5;
64  s = f * c * c;
65  for (i = 0; i < iMax; i++) {
66  t = (a + b) / 2;
67  c = (a - b) / 2;
68  b = sqrt (a * b);
69  a = t;
70  f *= 2;
71  s += f * c * c;
72  if (c / a < NR_EPSI) break;
73  }
74  if (i >= iMax) {
75  k = 0; e = 0;
76  }
77  else {
78  k = M_PI_2 / a;
79  e = M_PI_2 * (1 - s) / a;
80  if (arg < 0) {
81  k *= fk;
82  e *= fe;
83  }
84  }
85  }
86 }
87 
88 const nr_double_t SN_ACC = 1e-5; // Accuracy of sn(x) is SN_ACC^2
89 const nr_double_t K_ERR = 1e-8; // Accuracy of K(k)
90 
91 // Computes Carlson's elliptic integral of the first kind.
92 nr_double_t fspecial::ellip_rf (nr_double_t x, nr_double_t y, nr_double_t z) {
93  nr_double_t al, av, dx, dy, dz, e2, e3;
94  nr_double_t sx, sy, sz, xt, yt, zt;
95 
96  // constants
97  const nr_double_t c1 = 1.0 / 24.0;
98  const nr_double_t c2 = 0.1;
99  const nr_double_t c3 = 3.0 / 44.0;
100  const nr_double_t c4 = 1.0 / 14.0;
101 
102  xt = x;
103  yt = y;
104  zt = z;
105  do {
106  sx = sqrt (xt);
107  sy = sqrt (yt);
108  sz = sqrt (zt);
109  al = sx * (sy + sz) + sy * sz;
110  xt = 0.25 * (xt + al);
111  yt = 0.25 * (yt + al);
112  zt = 0.25 * (zt + al);
113  av = (xt + yt + zt) / 3.0;
114  dx = (av - xt) / av;
115  dy = (av - yt) / av;
116  dz = (av - zt) / av;
117  }
118  while (MAX (MAX (fabs (dx), fabs (dy)), fabs (dz)) > K_ERR);
119 
120  e2 = dx * dy - dz * dz;
121  e3 = dx * dy * dz;
122  return (1 + (c1 * e2 - c2 - c3 * e3) * e2 + c4 * e3) / sqrt (av);
123 }
124 
125 // Compute the Jacobian elliptic functions sn (u,k), cn (u,k) and dn (u,k).
126 nr_double_t fspecial::ellip_sncndn (nr_double_t u, nr_double_t k,
127  nr_double_t& sn, nr_double_t& cn, nr_double_t& dn) {
128  nr_double_t a, b, c, d;
129  nr_double_t fn[14], en[14];
130  int l;
131  bool bo;
132 
133  d = 1 - k;
134  if (k) {
135  bo = (k < 0);
136  if (bo) {
137  k /= -1 / d;
138  u *= (d = sqrt (d));
139  }
140  a = 1;
141  dn = 1;
142  for (int i = 1; i < 14; i++) {
143  l = i;
144  fn[i] = a;
145  en[i] = (k = sqrt (k));
146  c = 0.5 * (a + k);
147  if (fabs (a - k) <= SN_ACC * a)
148  break;
149  k *= a;
150  a = c;
151  }
152  u *= c;
153  sn = sin (u);
154  cn = cos (u);
155  if (sn) {
156  a = cn / sn;
157  c *= a;
158  for (int ii = l; ii > 0; ii--) {
159  b = fn[ii];
160  a *= c;
161  c *= dn;
162  dn = (en[ii] + a) / (b + a);
163  a = c / b;
164  }
165  a = 1 / sqrt (c * c + 1);
166  sn = (sn >= 0 ? a : -a);
167  cn = sn * c;
168  }
169  if (bo) {
170  a = dn;
171  dn = cn;
172  cn = a;
173  sn /= d;
174  }
175  }
176  else {
177  cn = 1 / cosh (u);
178  dn = cn;
179  sn = tanh (u);
180  }
181  return sn;
182 }
183 
184 /* data for a Chebyshev series over a given interval */
186  nr_double_t * c; /* coefficients */
187  int order; /* order of expansion */
188  nr_double_t a; /* lower interval point */
189  nr_double_t b; /* upper interval point */
190 };
191 typedef struct cheb_series_t cheb_series;
192 
193 static nr_double_t cheb_eval (const cheb_series * cs, const nr_double_t x) {
194  nr_double_t d = 0.0;
195  nr_double_t dd = 0.0;
196  nr_double_t y = (2.0 * x - cs->a - cs->b) / (cs->b - cs->a);
197  nr_double_t y2 = 2.0 * y;
198  for (int i = cs->order; i >= 1; i--) {
199  nr_double_t t = d;
200  d = y2 * d - dd + cs->c[i];
201  dd = t;
202  }
203  d = y * d - dd + 0.5 * cs->c[0];
204  return d;
205 }
206 
207 #if !defined (HAVE_ERF) || !defined (HAVE_ERFC)
208 
209 /* Chebyshev fit for erfc ((t+1)/2), -1 < t < 1 */
210 static nr_double_t erfc_xlt1_data[20] = {
211  1.06073416421769980345174155056,
212  -0.42582445804381043569204735291,
213  0.04955262679620434040357683080,
214  0.00449293488768382749558001242,
215  -0.00129194104658496953494224761,
216  -0.00001836389292149396270416979,
217  0.00002211114704099526291538556,
218  -5.23337485234257134673693179020e-7,
219  -2.78184788833537885382530989578e-7,
220  1.41158092748813114560316684249e-8,
221  2.72571296330561699984539141865e-9,
222  -2.06343904872070629406401492476e-10,
223  -2.14273991996785367924201401812e-11,
224  2.22990255539358204580285098119e-12,
225  1.36250074650698280575807934155e-13,
226  -1.95144010922293091898995913038e-14,
227  -6.85627169231704599442806370690e-16,
228  1.44506492869699938239521607493e-16,
229  2.45935306460536488037576200030e-18,
230  -9.29599561220523396007359328540e-19
231 };
233  erfc_xlt1_data, 19, -1, 1
234 };
235 
236 /* Chebyshev fit for erfc (x) * exp (x^2), 1 < x < 5, x = 2t + 3, -1 < t < 1 */
237 static nr_double_t erfc_x15_data[25] = {
238  0.44045832024338111077637466616,
239  -0.143958836762168335790826895326,
240  0.044786499817939267247056666937,
241  -0.013343124200271211203618353102,
242  0.003824682739750469767692372556,
243  -0.001058699227195126547306482530,
244  0.000283859419210073742736310108,
245  -0.000073906170662206760483959432,
246  0.000018725312521489179015872934,
247  -4.62530981164919445131297264430e-6,
248  1.11558657244432857487884006422e-6,
249  -2.63098662650834130067808832725e-7,
250  6.07462122724551777372119408710e-8,
251  -1.37460865539865444777251011793e-8,
252  3.05157051905475145520096717210e-9,
253  -6.65174789720310713757307724790e-10,
254  1.42483346273207784489792999706e-10,
255  -3.00141127395323902092018744545e-11,
256  6.22171792645348091472914001250e-12,
257  -1.26994639225668496876152836555e-12,
258  2.55385883033257575402681845385e-13,
259  -5.06258237507038698392265499770e-14,
260  9.89705409478327321641264227110e-15,
261  -1.90685978789192181051961024995e-15,
262  3.50826648032737849245113757340e-16
263 };
265  erfc_x15_data, 24, -1, 1
266 };
267 
268 /* Chebyshev fit for erfc(x) * exp(x^2),
269  5 < x < 10, x = (5t + 15)/2, -1 < t < */
270 static nr_double_t erfc_x510_data[20] = {
271  1.11684990123545698684297865808,
272  0.003736240359381998520654927536,
273  -0.000916623948045470238763619870,
274  0.000199094325044940833965078819,
275  -0.000040276384918650072591781859,
276  7.76515264697061049477127605790e-6,
277  -1.44464794206689070402099225301e-6,
278  2.61311930343463958393485241947e-7,
279  -4.61833026634844152345304095560e-8,
280  8.00253111512943601598732144340e-9,
281  -1.36291114862793031395712122089e-9,
282  2.28570483090160869607683087722e-10,
283  -3.78022521563251805044056974560e-11,
284  6.17253683874528285729910462130e-12,
285  -9.96019290955316888445830597430e-13,
286  1.58953143706980770269506726000e-13,
287  -2.51045971047162509999527428316e-14,
288  3.92607828989125810013581287560e-15,
289  -6.07970619384160374392535453420e-16,
290  9.12600607264794717315507477670e-17
291 };
293  erfc_x510_data, 19, -1, 1
294 };
295 
296 /* Estimates erfc (x) valid for 8 < x < 100, this is based on index 5725
297  in Hart et al. */
298 static nr_double_t erfc8 (nr_double_t x) {
299  static nr_double_t p[] = {
300  2.97886562639399288862,
301  7.409740605964741794425,
302  6.1602098531096305440906,
303  5.019049726784267463450058,
304  1.275366644729965952479585264,
305  0.5641895835477550741253201704
306  };
307  static nr_double_t q[] = {
308  3.3690752069827527677,
309  9.608965327192787870698,
310  17.08144074746600431571095,
311  12.0489519278551290360340491,
312  9.396034016235054150430579648,
313  2.260528520767326969591866945,
314  1.0
315  };
316  nr_double_t n = 0.0, d = 0.0;
317  int i;
318  n = p[5];
319  for (i = 4; i >= 0; --i) n = x * n + p[i];
320  d = q[6];
321  for (i = 5; i >= 0; --i) d = x * d + q[i];
322 
323  return exp (-x * x) * (n / d);
324 }
325 
326 #endif /* !HAVE_ERF || !HAVE_ERFC */
327 
328 #ifndef HAVE_ERFC
329 
330 nr_double_t fspecial::erfc (nr_double_t x) {
331  const nr_double_t ax = fabs (x);
332  nr_double_t val;
333 
334  if (ax <= 1.0) {
335  nr_double_t t = 2.0 * ax - 1.0;
336  val = cheb_eval (&erfc_xlt1_cs, t);
337  }
338  else if (ax <= 5.0) {
339  nr_double_t ex2 = exp (-x * x);
340  nr_double_t t = 0.5 * (ax - 3.0);
341  val = ex2 * cheb_eval (&erfc_x15_cs, t);
342  }
343  else if (ax < 10.0) {
344  nr_double_t ex = exp(-x * x) / ax;
345  nr_double_t t = (2.0 * ax - 15.0) / 5.0;
346  val = ex * cheb_eval (&erfc_x510_cs, t);
347  }
348  else {
349  val = erfc8 (ax);
350  }
351  return (x < 0.0) ? 2.0 - val : val;
352 }
353 #else
354 nr_double_t fspecial::erfc (nr_double_t d) {
355  return ::erfc (d);
356 }
357 #endif /* HAVE_ERFC */
358 
359 #ifndef HAVE_ERF
360 
361 /* Abramowitz + Stegun, 7.1.5 */
362 static nr_double_t erfseries (nr_double_t x) {
363  nr_double_t c = x;
364  nr_double_t e = c;
365  nr_double_t d;
366  for (int k = 1; k < 30; ++k) {
367  c *= -x * x / k;
368  d = c / (2.0 * k + 1.0);
369  e += d;
370  }
371  return 2.0 / M_SQRTPI * e;
372 }
373 
374 nr_double_t fspecial::erf (nr_double_t x) {
375  if (fabs (x) < 1.0) {
376  return erfseries (x);
377  }
378  return 1.0 - erfc (x);
379 }
380 
381 #else
382 nr_double_t fspecial::erf (nr_double_t d) {
383  return ::erf (d);
384 }
385 #endif /* HAVE_ERF */
386 
387 // Inverse of the error function erf().
388 nr_double_t fspecial::erfinv (nr_double_t y) {
389  nr_double_t x = 0.0; // The output
390  nr_double_t z = 0.0; // Intermediate variable
391  nr_double_t y0 = 0.7; // Central range variable
392 
393  // Coefficients in rational approximations.
394  nr_double_t a[4] = { 0.886226899, -1.645349621, 0.914624893, -0.140543331};
395  nr_double_t b[4] = {-2.118377725, 1.442710462, -0.329097515, 0.012229801};
396  nr_double_t c[4] = {-1.970840454, -1.624906493, 3.429567803, 1.641345311};
397  nr_double_t d[2] = { 3.543889200, 1.637067800};
398 
399  if (y < -1.0 || 1.0 < y) {
400  x = log (-1.0);
401  }
402  else if (y == -1.0 || 1.0 == y) {
403  x = -y * log(0.0);
404  }
405  else if (-1.0 < y && y < -y0) {
406  z = sqrt(-log((1.0+y)/2.0));
407  x = -(((c[3]*z+c[2])*z+c[1])*z+c[0])/((d[1]*z+d[0])*z+1.0);
408  }
409  else {
410  if (-y0 < y && y < y0) {
411  z = y * y;
412  x = y*(((a[3]*z+a[2])*z+a[1])*z+a[0]) /
413  ((((b[3]*z+b[3])*z+b[1])*z+b[0])*z+1.0);
414  }
415  else if (y0 < y && y < 1.0) {
416  z = sqrt(-log((1.0-y)/2.0));
417  x = (((c[3]*z+c[2])*z+c[1])*z+c[0])/((d[1]*z+d[0])*z+1.0);
418  }
419 
420  // Two steps of Newton-Raphson correction to full accuracy.
421  x = x - (erf (x) - y) / (2.0 / M_SQRTPI * exp (-x * x));
422  x = x - (erf (x) - y) / (2.0 / M_SQRTPI * exp (-x * x));
423  }
424  return x;
425 }
426 
427 static nr_double_t bi0_data[12] = {
428  -.07660547252839144951,
429  1.92733795399380827000,
430  .22826445869203013390,
431  .01304891466707290428,
432  .00043442709008164874,
433  .00000942265768600193,
434  .00000014340062895106,
435  .00000000161384906966,
436  .00000000001396650044,
437  .00000000000009579451,
438  .00000000000000053339,
439  .00000000000000000245
440 };
441 static cheb_series bi0_cs = {
442  bi0_data, 11, -1, 1
443 };
444 
445 static nr_double_t ai0_data[21] = {
446  .07575994494023796,
447  .00759138081082334,
448  .00041531313389237,
449  .00001070076463439,
450  -.00000790117997921,
451  -.00000078261435014,
452  .00000027838499429,
453  .00000000825247260,
454  -.00000001204463945,
455  .00000000155964859,
456  .00000000022925563,
457  -.00000000011916228,
458  .00000000001757854,
459  .00000000000112822,
460  -.00000000000114684,
461  .00000000000027155,
462  -.00000000000002415,
463  -.00000000000000608,
464  .00000000000000314,
465  -.00000000000000071,
466  .00000000000000007
467 };
468 static cheb_series ai0_cs = {
469  ai0_data, 20, -1, 1
470 };
471 
472 static nr_double_t ai02_data[22] = {
473  .05449041101410882,
474  .00336911647825569,
475  .00006889758346918,
476  .00000289137052082,
477  .00000020489185893,
478  .00000002266668991,
479  .00000000339623203,
480  .00000000049406022,
481  .00000000001188914,
482  -.00000000003149915,
483  -.00000000001321580,
484  -.00000000000179419,
485  .00000000000071801,
486  .00000000000038529,
487  .00000000000001539,
488  -.00000000000004151,
489  -.00000000000000954,
490  .00000000000000382,
491  .00000000000000176,
492  -.00000000000000034,
493  -.00000000000000027,
494  .00000000000000003
495 };
497  ai02_data, 21, -1, 1
498 };
499 
500 // Modified Bessel function of order zero.
501 nr_double_t fspecial::i0 (nr_double_t x) {
502  nr_double_t y = fabs (x);
503  nr_double_t val;
504 
505  if (y < 2.0 * sqrt (NR_EPSI)) {
506  val = 1.0;
507  }
508  else if (y <= 3.0) {
509  val = 2.75 + cheb_eval (&bi0_cs, y * y / 4.5 - 1.0);
510  }
511  else if (y <= 8.0) {
512  val = cheb_eval (&ai0_cs, (48.0 / y - 11.0) / 5.0);
513  val = exp (y) * (0.375 + val) / sqrt (y);
514  }
515  else {
516  val = cheb_eval (&ai02_cs, 16.0 / y - 1.0);
517  val = exp (y) * (0.375 + val) / sqrt (y);
518  }
519  return val;
520 }
521 
522 // Lower tail quantile for the standard normal distribution function.
523 nr_double_t fspecial::ltqnorm (nr_double_t x) {
524  nr_double_t q, r, e, u, z = 0.0;
525  static nr_double_t a[] = {
526  -3.969683028665376e+01, 2.209460984245205e+02,
527  -2.759285104469687e+02, 1.383577518672690e+02,
528  -3.066479806614716e+01, 2.506628277459239e+00 };
529  static nr_double_t b[] = {
530  -5.447609879822406e+01, 1.615858368580409e+02,
531  -1.556989798598866e+02, 6.680131188771972e+01,
532  -1.328068155288572e+01 };
533  static nr_double_t c[] = {
534  -7.784894002430293e-03, -3.223964580411365e-01,
535  -2.400758277161838e+00, -2.549732539343734e+00,
536  4.374664141464968e+00, 2.938163982698783e+00 };
537  static nr_double_t d[] = {
538  7.784695709041462e-03, 3.224671290700398e-01,
539  2.445134137142996e+00, 3.754408661907416e+00 };
540 
541  // Define break-points.
542  nr_double_t pl = 0.02425;
543  nr_double_t ph = 1.0 - pl;
544 
545  // Rational approximation for central region:
546  if (pl <= x && x <= ph) {
547  q = x - 0.5;
548  r = q * q;
549  z = (((((a[0]*r+a[1])*r+a[2])*r+a[3])*r+a[4])*r+a[5])*q/
550  (((((b[0]*r+b[1])*r+b[2])*r+b[3])*r+b[4])*r+1);
551  }
552  // Rational approximation for lower region:
553  else if (0.0 < x && x < pl) {
554  q = sqrt(-2*log(x));
555  z = (((((c[0]*q+c[1])*q+c[2])*q+c[3])*q+c[4])*q+c[5])/
556  ((((d[0]*q+d[1])*q+d[2])*q+d[3])*q+1);
557  }
558  // Rational approximation for upper region:
559  else if (ph < x && x < 1.0) {
560  q = sqrt(-2*log(1-x));
561  z = -(((((c[0]*q+c[1])*q+c[2])*q+c[3])*q+c[4])*q+c[5])/
562  ((((d[0]*q+d[1])*q+d[2])*q+d[3])*q+1);
563  }
564  // Case when X = 0:
565  else if (x == 0.0) {
566  z = -NR_INF;
567  }
568  // Case when X = 1:
569  else if (x == 1.0) {
570  z = +NR_INF;
571  }
572  // Cases when output will be NaN:
573  else if (x < 0.0 || x > 1.0 || std::isnan (x)) {
574  z = +NR_NAN;
575  }
576 
577  // The relative error of the approximation has absolute value less
578  // than 1.15e-9. One iteration of Halley's rational method (third
579  // order) gives full machine precision.
580  if (0.0 < x && x < 1.0) {
581  e = 0.5 * erfc (-z / M_SQRT2) - x; // error
582  u = e * M_SQRT2 * M_SQRTPI * exp (z * z / 2); // f(z)/df(z)
583  z = z - u / (1 + z * u / 2); // Halley's method
584  }
585  return z;
586 }
587 
588 // Inverse of the error function erfc().
589 nr_double_t fspecial::erfcinv (nr_double_t x) {
590  return -ltqnorm (x / 2.0) / M_SQRT2;
591 }
592 
static nr_double_t erfc_xlt1_data[20]
Definition: fspecial.cpp:210
nr_double_t erf(nr_double_t)
Definition: fspecial.cpp:374
l
Definition: parse_vcd.y:213
nr_double_t ltqnorm(nr_double_t)
Definition: fspecial.cpp:523
nr_double_t ellip_sncndn(nr_double_t, nr_double_t, nr_double_t &, nr_double_t &, nr_double_t &)
Definition: fspecial.cpp:126
nr_complex_t erf(const nr_complex_t z)
Error function.
Definition: complex.cpp:766
nr_double_t * c
Definition: fspecial.cpp:186
const nr_double_t K_ERR
Definition: fspecial.cpp:89
nr_double_t erfinv(nr_double_t)
Definition: fspecial.cpp:388
static nr_double_t erfseries(nr_double_t x)
Definition: fspecial.cpp:362
static cheb_series ai02_cs
Definition: fspecial.cpp:496
nr_complex_t cos(const nr_complex_t z)
Compute complex cosine.
Definition: complex.cpp:57
t
Definition: parse_vcd.y:290
static cheb_series bi0_cs
Definition: fspecial.cpp:441
Global physical constants header file.
nr_double_t erfcinv(nr_double_t)
Definition: fspecial.cpp:589
nr_double_t b
Definition: fspecial.cpp:189
n
Definition: parse_citi.y:147
r
Definition: parse_mdl.y:515
static cheb_series erfc_x15_cs
Definition: fspecial.cpp:264
nr_complex_t cosh(const nr_complex_t z)
Compute complex hyperbolic cosine.
Definition: complex.cpp:135
i
Definition: parse_mdl.y:516
#define M_PI_2
Half of Archimedes' constant ( )
Definition: consts.h:51
#define NR_INF
Definition: precision.h:113
nr_complex_t sqrt(const nr_complex_t z)
Compute principal value of square root.
Definition: complex.cpp:271
static nr_double_t bi0_data[12]
Definition: fspecial.cpp:427
const nr_double_t SN_ACC
Definition: fspecial.cpp:88
static nr_double_t cheb_eval(const cheb_series *cs, const nr_double_t x)
Definition: fspecial.cpp:193
nr_complex_t tanh(const nr_complex_t z)
Compute complex hyperbolic tangent.
Definition: complex.cpp:153
nr_complex_t sin(const nr_complex_t z)
Compute complex sine.
Definition: complex.cpp:66
eqn::constant * c1
static nr_double_t ai0_data[21]
Definition: fspecial.cpp:445
eqn::constant * c2
x
Definition: parse_mdl.y:498
void ellip_ke(nr_double_t, nr_double_t &, nr_double_t &)
Definition: fspecial.cpp:42
static nr_double_t ai02_data[22]
Definition: fspecial.cpp:472
#define M_SQRTPI
Square root of Archimedes' constant ( )
Definition: consts.h:67
nr_complex_t erfc(const nr_complex_t z)
Complementart error function.
Definition: complex.cpp:784
static cheb_series erfc_x510_cs
Definition: fspecial.cpp:292
nr_double_t a
Definition: fspecial.cpp:188
y
Definition: parse_mdl.y:499
static nr_double_t erfc_x15_data[25]
Definition: fspecial.cpp:237
nr_double_t erfc(nr_double_t)
Definition: fspecial.cpp:330
nr_complex_t exp(const nr_complex_t z)
Compute complex exponential.
Definition: complex.cpp:205
static cheb_series ai0_cs
Definition: fspecial.cpp:468
nr_complex_t log(const nr_complex_t z)
Compute principal value of natural logarithm of z.
Definition: complex.cpp:215
#define M_SQRT2
Square root of 2 ( )
Definition: consts.h:91
static nr_double_t erfc_x510_data[20]
Definition: fspecial.cpp:270
nr_double_t ellip_rf(nr_double_t, nr_double_t, nr_double_t)
Definition: fspecial.cpp:92
#define NR_NAN
Definition: precision.h:114
static nr_double_t erfc8(nr_double_t x)
Definition: fspecial.cpp:298
nr_double_t i0(nr_double_t)
Definition: fspecial.cpp:501
matrix arg(matrix a)
Computes the argument of each matrix element.
Definition: matrix.cpp:555
static cheb_series erfc_xlt1_cs
Definition: fspecial.cpp:232
#define MAX(x, y)
Maximum of x and y.
Definition: constants.h:127