source: sasmodels/sasmodels/special.py @ 0a9fcab

Last change on this file since 0a9fcab was 0a9fcab, checked in by Paul Kienzle <pkienzle@…>, 7 years ago

simplify handling of gauss.z in py2c

  • Property mode set to 100644
File size: 23.6 KB
Line 
1r"""
2Special Functions
3.................
4
5This following standard C99 math functions are available:
6
7    M_PI, M_PI_2, M_PI_4, M_SQRT1_2, M_E:
8        $\pi$, $\pi/2$, $\pi/4$, $1/\sqrt{2}$ and Euler's constant $e$
9
10    exp, log, pow(x,y), expm1, log1p, sqrt, cbrt:
11        Power functions $e^x$, $\ln x$, $x^y$, $e^x - 1$, $\ln 1 + x$,
12        $\sqrt{x}$, $\sqrt[3]{x}$. The functions expm1(x) and log1p(x)
13        are accurate across all $x$, including $x$ very close to zero.
14
15    sin, cos, tan, asin, acos, atan:
16        Trigonometry functions and inverses, operating on radians.
17
18    sinh, cosh, tanh, asinh, acosh, atanh:
19        Hyperbolic trigonometry functions.
20
21    atan2(y,x):
22        Angle from the $x$\ -axis to the point $(x,y)$, which is equal to
23        $\tan^{-1}(y/x)$ corrected for quadrant.  That is, if $x$ and $y$ are
24        both negative, then atan2(y,x) returns a value in quadrant III where
25        atan(y/x) would return a value in quadrant I. Similarly for
26        quadrants II and IV when $x$ and $y$ have opposite sign.
27
28    fabs(x), fmin(x,y), fmax(x,y), trunc, rint:
29        Floating point functions.  rint(x) returns the nearest integer.
30
31    NAN:
32        NaN, Not a Number, $0/0$.  Use isnan(x) to test for NaN.  Note that
33        you cannot use :code:`x == NAN` to test for NaN values since that
34        will always return false.  NAN does not equal NAN!  The alternative,
35        :code:`x != x` may fail if the compiler optimizes the test away.
36
37    INFINITY:
38        $\infty, 1/0$.  Use isinf(x) to test for infinity, or isfinite(x)
39        to test for finite and not NaN.
40
41    erf, erfc, tgamma, lgamma:  **do not use**
42        Special functions that should be part of the standard, but are missing
43        or inaccurate on some platforms. Use sas_erf, sas_erfc and sas_gamma
44        instead (see below). Note: lgamma(x) has not yet been tested.
45
46Some non-standard constants and functions are also provided:
47
48    M_PI_180, M_4PI_3:
49        $\frac{\pi}{180}$, $\frac{4\pi}{3}$
50
51    SINCOS(x, s, c):
52        Macro which sets s=sin(x) and c=cos(x). The variables *c* and *s*
53        must be declared first.
54
55    square(x):
56        $x^2$
57
58    cube(x):
59        $x^3$
60
61    sas_sinx_x(x):
62        $\sin(x)/x$, with limit $\sin(0)/0 = 1$.
63
64    powr(x, y):
65        $x^y$ for $x \ge 0$; this is faster than general $x^y$ on some GPUs.
66
67    pown(x, n):
68        $x^n$ for $n$ integer; this is faster than general $x^n$ on some GPUs.
69
70    FLOAT_SIZE:
71        The number of bytes in a floating point value.  Even though all
72        variables are declared double, they may be converted to single
73        precision float before running. If your algorithm depends on
74        precision (which is not uncommon for numerical algorithms), use
75        the following::
76
77            #if FLOAT_SIZE>4
78            ... code for double precision ...
79            #else
80            ... code for single precision ...
81            #endif
82
83    SAS_DOUBLE:
84        A replacement for :code:`double` so that the declared variable will
85        stay double precision; this should generally not be used since some
86        graphics cards do not support double precision.  There is no provision
87        for forcing a constant to stay double precision.
88
89The following special functions and scattering calculations are defined.
90These functions have been tuned to be fast and numerically stable down
91to $q=0$ even in single precision.  In some cases they work around bugs
92which appear on some platforms but not others, so use them where needed.
93Add the files listed in :code:`source = ["lib/file.c", ...]` to your *model.py*
94file in the order given, otherwise these functions will not be available.
95
96    polevl(x, c, n):
97        Polynomial evaluation $p(x) = \sum_{i=0}^n c_i x^i$ using Horner's
98        method so it is faster and more accurate.
99
100        $c = \{c_n, c_{n-1}, \ldots, c_0 \}$ is the table of coefficients,
101        sorted from highest to lowest.
102
103    p1evl(x, c, n):
104        Evaluate normalized polynomial $p(x) = x^n + \sum_{i=0}^{n-1} c_i x^i$
105        using Horner's method so it is faster and more accurate.
106
107        $c = \{c_{n-1}, c_{n-2} \ldots, c_0 \}$ is the table of coefficients,
108        sorted from highest to lowest.
109
110    sas_gamma(x):
111        Gamma function $\text{sas_gamma}(x) = \Gamma(x)$.
112
113        The standard math function, tgamma(x) is unstable for $x < 1$
114        on some platforms.
115
116    sas_erf(x), sas_erfc(x):
117        Error function
118        $\text{sas_erf}(x) = \frac{2}{\sqrt\pi}\int_0^x e^{-t^2}\,dt$
119        and complementary error function
120        $\text{sas_erfc}(x) = \frac{2}{\sqrt\pi}\int_x^{\infty} e^{-t^2}\,dt$.
121
122        The standard math functions erf(x) and erfc(x) are slower and broken
123        on some platforms.
124
125    sas_J0(x):
126        Bessel function of the first kind $\text{sas_J0}(x)=J_0(x)$ where
127        $J_0(x) = \frac{1}{\pi}\int_0^\pi \cos(x\sin(\tau))\,d\tau$.
128
129        The standard math function j0(x) is not available on all platforms.
130
131    sas_J1(x):
132        Bessel function of the first kind  $\text{sas_J1}(x)=J_1(x)$ where
133        $J_1(x) = \frac{1}{\pi}\int_0^\pi \cos(\tau - x\sin(\tau))\,d\tau$.
134
135        The standard math function j1(x) is not available on all platforms.
136
137    sas_JN(n, x):
138        Bessel function of the first kind and integer order $n$:
139        $\text{sas_JN}(n, x)=J_n(x)$ where
140        $J_n(x) = \frac{1}{\pi}\int_0^\pi \cos(n\tau - x\sin(\tau))\,d\tau$.
141        If $n$ = 0 or 1, it uses sas_J0(x) or sas_J1(x), respectively.
142
143        The standard math function jn(n, x) is not available on all platforms.
144
145    sas_Si(x):
146        Sine integral $\text{Si}(x) = \int_0^x \tfrac{\sin t}{t}\,dt$.
147
148        This function uses Taylor series for small and large arguments:
149
150        For large arguments,
151
152        .. math::
153
154            \text{Si}(x) \sim \frac{\pi}{2}
155             - \frac{\cos(x)}{x}
156            \left(1 - \frac{2!}{x^2} + \frac{4!}{x^4} - \frac{6!}{x^6} \right)
157             - \frac{\sin(x)}{x}
158            \left(\frac{1}{x} - \frac{3!}{x^3} + \frac{5!}{x^5} - \frac{7!}{x^7}\right)
159
160        For small arguments,
161
162        .. math::
163
164           \text{Si}(x) \sim x
165           - \frac{x^3}{3\times 3!} + \frac{x^5}{5 \times 5!} - \frac{x^7}{7 \times 7!}
166           + \frac{x^9}{9\times 9!} - \frac{x^{11}}{11\times 11!}
167
168    sas_3j1x_x(x):
169        Spherical Bessel form
170        $\text{sph_j1c}(x) = 3 j_1(x)/x = 3 (\sin(x) - x \cos(x))/x^3$,
171        with a limiting value of 1 at $x=0$, where $j_1(x)$ is the spherical
172        Bessel function of the first kind and first order.
173
174        This function uses a Taylor series for small $x$ for numerical accuracy.
175
176
177    sas_2J1x_x(x):
178        Bessel form $\text{sas_J1c}(x) = 2 J_1(x)/x$, with a limiting value
179        of 1 at $x=0$, where $J_1(x)$ is the Bessel function of first kind
180        and first order.
181
182
183    gauss76.n, gauss76.z[i], gauss76.w[i]:
184        Points $z_i$ and weights $w_i$ for 76-point Gaussian quadrature, respectively,
185        computing $\int_{-1}^1 f(z)\,dz \approx \sum_{i=1}^{76} w_i\,f(z_i)$.
186        When translating the model to C, include 'lib/gauss76.c' in the source
187        and use :code:`GAUSS_N`, :code:`GAUSS_Z`, and :code:`GAUSS_W`.
188
189        Similar arrays are available in :code:`gauss20` for 20-point quadrature
190        and :code:`gauss150.c` for 150-point quadrature. By using
191        :code:`import gauss76 as gauss` it is easy to change the number of
192        points in the integration.
193"""
194# pylint: disable=unused-import
195
196import numpy as np
197
198# Functions to add to our standard set
199from numpy import degrees, radians
200
201# C99 standard math library functions
202from numpy import exp, log, power as pow, expm1, logp1, sqrt, cbrt
203from numpy import sin, cos, tan, arcsin as asin, arccos as acos, arctan as atan
204from numpy import sinh, cosh, tanh, arcsinh as asinh, arccosh as acosh, arctanh as atanh
205from numpy import arctan2 as atan2
206from numpy import fabs, fmin, fmax, trunc, rint
207from numpy import pi, nan, inf
208from scipy.special import gamma as sas_gamma
209from scipy.special import erf as sas_erf
210from scipy.special import erfc as sas_erfc
211from scipy.special import j0 as sas_J0
212from scipy.special import j1 as sas_J1
213from scipy.special import jn as sas_JN
214
215# erf, erfc, tgamma, lgamma  **do not use**
216
217# C99 standard math constants
218M_PI, M_PI_2, M_PI_4, M_SQRT1_2, M_E = np.pi, np.pi/2, np.pi/4, np.sqrt(0.5), np.e
219NAN = nan
220INFINITY = inf
221
222# non-standard constants
223M_PI_180, M_4PI_3 = M_PI/180, 4*M_PI/3
224
225# can't do SINCOS in python; use "s, c = SINCOS(x)" instead
226def SINCOS(x):
227    """return sin(x), cos(x)"""
228    return sin(x), cos(x)
229sincos = SINCOS
230
231def square(x):
232    """return x^2"""
233    return x*x
234
235def cube(x):
236    """return x^3"""
237    return x*x*x
238
239def sas_sinx_x(x):
240    """return sin(x)/x"""
241    from numpy import sinc as _sinc
242    return _sinc(x/M_PI)
243
244def powr(x, y):
245    """return x^y for x>0"""
246    return x**y
247def pown(x, n):
248    """return x^n for n integer"""
249    return x**n
250
251FLOAT_SIZE = 8
252
253def polevl(x, c, n):
254    """return p(x) for polynomial p of degree n-1 with coefficients c"""
255    return np.polyval(c[:n], x)
256
257def p1evl(x, c, n):
258    """return x^n + p(x) for polynomial p of degree n-1 with coefficients c"""
259    return np.polyval(np.hstack(([1.], c))[:n], x)
260
261def sas_Si(x):
262    """return Si(x)"""
263    from scipy.special import sici
264    return sici(x)[0]
265
266def sas_j1(x):
267    """return j1(x)"""
268    if np.isscalar(x):
269        retvalue = (sin(x) - x*cos(x))/x**2 if x != 0. else 0.
270    else:
271        with np.errstate(all='ignore'):
272            retvalue = (sin(x) - x*cos(x))/x**2
273        retvalue[x == 0.] = 0.
274    return retvalue
275
276def sas_3j1x_x(x):
277    """return 3*j1(x)/x"""
278    if np.isscalar(x):
279        retvalue = 3*(sin(x) - x*cos(x))/x**3 if x != 0. else 1.
280    else:
281        with np.errstate(all='ignore'):
282            retvalue = 3*(sin(x) - x*cos(x))/x**3
283        retvalue[x == 0.] = 1.
284    return retvalue
285
286def sas_2J1x_x(x):
287    """return 2*J1(x)/x"""
288    if np.isscalar(x):
289        retvalue = 2*sas_J1(x)/x if x != 0 else 1.
290    else:
291        with np.errstate(all='ignore'):
292            retvalue = 2*sas_J1(x)/x
293        retvalue[x == 0] = 1.
294    return retvalue
295
296
297# Gaussians
298class Gauss:
299    def __init__(self, w, z):
300        self.n = len(w)
301        self.w = w
302        self.z = z
303
304gauss20 = Gauss(
305    w=np.array([
306        .0176140071391521,
307        .0406014298003869,
308        .0626720483341091,
309        .0832767415767047,
310        .10193011981724,
311        .118194531961518,
312        .131688638449177,
313        .142096109318382,
314        .149172986472604,
315        .152753387130726,
316        .152753387130726,
317        .149172986472604,
318        .142096109318382,
319        .131688638449177,
320        .118194531961518,
321        .10193011981724,
322        .0832767415767047,
323        .0626720483341091,
324        .0406014298003869,
325        .0176140071391521
326    ]),
327    z=np.array([
328        -.993128599185095,
329        -.963971927277914,
330        -.912234428251326,
331        -.839116971822219,
332        -.746331906460151,
333        -.636053680726515,
334        -.510867001950827,
335        -.37370608871542,
336        -.227785851141645,
337        -.076526521133497,
338        .0765265211334973,
339        .227785851141645,
340        .37370608871542,
341        .510867001950827,
342        .636053680726515,
343        .746331906460151,
344        .839116971822219,
345        .912234428251326,
346        .963971927277914,
347        .993128599185095
348    ])
349)
350
351gauss76 = Gauss(
352    w=np.array([
353        .00126779163408536,             #0
354        .00294910295364247,
355        .00462793522803742,
356        .00629918049732845,
357        .00795984747723973,
358        .00960710541471375,
359        .0112381685696677,
360        .0128502838475101,
361        .0144407317482767,
362        .0160068299122486,
363        .0175459372914742,              #10
364        .0190554584671906,
365        .020532847967908,
366        .0219756145344162,
367        .0233813253070112,
368        .0247476099206597,
369        .026072164497986,
370        .0273527555318275,
371        .028587223650054,
372        .029773487255905,
373        .0309095460374916,              #20
374        .0319934843404216,
375        .0330234743977917,
376        .0339977794120564,
377        .0349147564835508,
378        .0357728593807139,
379        .0365706411473296,
380        .0373067565423816,
381        .0379799643084053,
382        .0385891292645067,
383        .0391332242205184,              #30
384        .0396113317090621,
385        .0400226455325968,
386        .040366472122844,
387        .0406422317102947,
388        .0408494593018285,
389        .040987805464794,
390        .0410570369162294,
391        .0410570369162294,
392        .040987805464794,
393        .0408494593018285,              #40
394        .0406422317102947,
395        .040366472122844,
396        .0400226455325968,
397        .0396113317090621,
398        .0391332242205184,
399        .0385891292645067,
400        .0379799643084053,
401        .0373067565423816,
402        .0365706411473296,
403        .0357728593807139,              #50
404        .0349147564835508,
405        .0339977794120564,
406        .0330234743977917,
407        .0319934843404216,
408        .0309095460374916,
409        .029773487255905,
410        .028587223650054,
411        .0273527555318275,
412        .026072164497986,
413        .0247476099206597,              #60
414        .0233813253070112,
415        .0219756145344162,
416        .020532847967908,
417        .0190554584671906,
418        .0175459372914742,
419        .0160068299122486,
420        .0144407317482767,
421        .0128502838475101,
422        .0112381685696677,
423        .00960710541471375,             #70
424        .00795984747723973,
425        .00629918049732845,
426        .00462793522803742,
427        .00294910295364247,
428        .00126779163408536              #75 (indexed from 0)
429    ]),
430    z=np.array([
431        -.999505948362153,              #0
432        -.997397786355355,
433        -.993608772723527,
434        -.988144453359837,
435        -.981013938975656,
436        -.972229228520377,
437        -.961805126758768,
438        -.949759207710896,
439        -.936111781934811,
440        -.92088586125215,
441        -.904107119545567,              #10
442        -.885803849292083,
443        -.866006913771982,
444        -.844749694983342,
445        -.822068037328975,
446        -.7980001871612,
447        -.77258672828181,
448        -.74587051350361,
449        -.717896592387704,
450        -.688712135277641,
451        -.658366353758143,              #20
452        -.626910417672267,
453        -.594397368836793,
454        -.560882031601237,
455        -.526420920401243,
456        -.491072144462194,
457        -.454895309813726,
458        -.417951418780327,
459        -.380302767117504,
460        -.342012838966962,
461        -.303146199807908,              #30
462        -.263768387584994,
463        -.223945802196474,
464        -.183745593528914,
465        -.143235548227268,
466        -.102483975391227,
467        -.0615595913906112,
468        -.0205314039939986,
469        .0205314039939986,
470        .0615595913906112,
471        .102483975391227,                       #40
472        .143235548227268,
473        .183745593528914,
474        .223945802196474,
475        .263768387584994,
476        .303146199807908,
477        .342012838966962,
478        .380302767117504,
479        .417951418780327,
480        .454895309813726,
481        .491072144462194,               #50
482        .526420920401243,
483        .560882031601237,
484        .594397368836793,
485        .626910417672267,
486        .658366353758143,
487        .688712135277641,
488        .717896592387704,
489        .74587051350361,
490        .77258672828181,
491        .7980001871612, #60
492        .822068037328975,
493        .844749694983342,
494        .866006913771982,
495        .885803849292083,
496        .904107119545567,
497        .92088586125215,
498        .936111781934811,
499        .949759207710896,
500        .961805126758768,
501        .972229228520377,               #70
502        .981013938975656,
503        .988144453359837,
504        .993608772723527,
505        .997397786355355,
506        .999505948362153                #75
507    ])
508)
509
510gauss150 = Gauss(
511    z=np.array([
512        -0.9998723404457334,
513        -0.9993274305065947,
514        -0.9983473449340834,
515        -0.9969322929775997,
516        -0.9950828645255290,
517        -0.9927998590434373,
518        -0.9900842691660192,
519        -0.9869372772712794,
520        -0.9833602541697529,
521        -0.9793547582425894,
522        -0.9749225346595943,
523        -0.9700655145738374,
524        -0.9647858142586956,
525        -0.9590857341746905,
526        -0.9529677579610971,
527        -0.9464345513503147,
528        -0.9394889610042837,
529        -0.9321340132728527,
530        -0.9243729128743136,
531        -0.9162090414984952,
532        -0.9076459563329236,
533        -0.8986873885126239,
534        -0.8893372414942055,
535        -0.8795995893549102,
536        -0.8694786750173527,
537        -0.8589789084007133,
538        -0.8481048644991847,
539        -0.8368612813885015,
540        -0.8252530581614230,
541        -0.8132852527930605,
542        -0.8009630799369827,
543        -0.7882919086530552,
544        -0.7752772600680049,
545        -0.7619248049697269,
546        -0.7482403613363824,
547        -0.7342298918013638,
548        -0.7198995010552305,
549        -0.7052554331857488,
550        -0.6903040689571928,
551        -0.6750519230300931,
552        -0.6595056411226444,
553        -0.6436719971150083,
554        -0.6275578900977726,
555        -0.6111703413658551,
556        -0.5945164913591590,
557        -0.5776035965513142,
558        -0.5604390262878617,
559        -0.5430302595752546,
560        -0.5253848818220803,
561        -0.5075105815339176,
562        -0.4894151469632753,
563        -0.4711064627160663,
564        -0.4525925063160997,
565        -0.4338813447290861,
566        -0.4149811308476706,
567        -0.3959000999390257,
568        -0.3766465660565522,
569        -0.3572289184172501,
570        -0.3376556177463400,
571        -0.3179351925907259,
572        -0.2980762356029071,
573        -0.2780873997969574,
574        -0.2579773947782034,
575        -0.2377549829482451,
576        -0.2174289756869712,
577        -0.1970082295132342,
578        -0.1765016422258567,
579        -0.1559181490266516,
580        -0.1352667186271445,
581        -0.1145563493406956,
582        -0.0937960651617229,
583        -0.0729949118337358,
584        -0.0521619529078925,
585        -0.0313062657937972,
586        -0.0104369378042598,
587        0.0104369378042598,
588        0.0313062657937972,
589        0.0521619529078925,
590        0.0729949118337358,
591        0.0937960651617229,
592        0.1145563493406956,
593        0.1352667186271445,
594        0.1559181490266516,
595        0.1765016422258567,
596        0.1970082295132342,
597        0.2174289756869712,
598        0.2377549829482451,
599        0.2579773947782034,
600        0.2780873997969574,
601        0.2980762356029071,
602        0.3179351925907259,
603        0.3376556177463400,
604        0.3572289184172501,
605        0.3766465660565522,
606        0.3959000999390257,
607        0.4149811308476706,
608        0.4338813447290861,
609        0.4525925063160997,
610        0.4711064627160663,
611        0.4894151469632753,
612        0.5075105815339176,
613        0.5253848818220803,
614        0.5430302595752546,
615        0.5604390262878617,
616        0.5776035965513142,
617        0.5945164913591590,
618        0.6111703413658551,
619        0.6275578900977726,
620        0.6436719971150083,
621        0.6595056411226444,
622        0.6750519230300931,
623        0.6903040689571928,
624        0.7052554331857488,
625        0.7198995010552305,
626        0.7342298918013638,
627        0.7482403613363824,
628        0.7619248049697269,
629        0.7752772600680049,
630        0.7882919086530552,
631        0.8009630799369827,
632        0.8132852527930605,
633        0.8252530581614230,
634        0.8368612813885015,
635        0.8481048644991847,
636        0.8589789084007133,
637        0.8694786750173527,
638        0.8795995893549102,
639        0.8893372414942055,
640        0.8986873885126239,
641        0.9076459563329236,
642        0.9162090414984952,
643        0.9243729128743136,
644        0.9321340132728527,
645        0.9394889610042837,
646        0.9464345513503147,
647        0.9529677579610971,
648        0.9590857341746905,
649        0.9647858142586956,
650        0.9700655145738374,
651        0.9749225346595943,
652        0.9793547582425894,
653        0.9833602541697529,
654        0.9869372772712794,
655        0.9900842691660192,
656        0.9927998590434373,
657        0.9950828645255290,
658        0.9969322929775997,
659        0.9983473449340834,
660        0.9993274305065947,
661        0.9998723404457334
662    ]),
663    w=np.array([
664        0.0003276086705538,
665        0.0007624720924706,
666        0.0011976474864367,
667        0.0016323569986067,
668        0.0020663664924131,
669        0.0024994789888943,
670        0.0029315036836558,
671        0.0033622516236779,
672        0.0037915348363451,
673        0.0042191661429919,
674        0.0046449591497966,
675        0.0050687282939456,
676        0.0054902889094487,
677        0.0059094573005900,
678        0.0063260508184704,
679        0.0067398879387430,
680        0.0071507883396855,
681        0.0075585729801782,
682        0.0079630641773633,
683        0.0083640856838475,
684        0.0087614627643580,
685        0.0091550222717888,
686        0.0095445927225849,
687        0.0099300043714212,
688        0.0103110892851360,
689        0.0106876814158841,
690        0.0110596166734735,
691        0.0114267329968529,
692        0.0117888704247183,
693        0.0121458711652067,
694        0.0124975796646449,
695        0.0128438426753249,
696        0.0131845093222756,
697        0.0135194311690004,
698        0.0138484622795371,
699        0.0141714592928592,
700        0.0144882814685445,
701        0.0147987907597169,
702        0.0151028518701744,
703        0.0154003323133401,
704        0.0156911024699895,
705        0.0159750356447283,
706        0.0162520081211971,
707        0.0165218992159766,
708        0.0167845913311726,
709        0.0170399700056559,
710        0.0172879239649355,
711        0.0175283451696437,
712        0.0177611288626114,
713        0.0179861736145128,
714        0.0182033813680609,
715        0.0184126574807331,
716        0.0186139107660094,
717        0.0188070535331042,
718        0.0189920016251754,
719        0.0191686744559934,
720        0.0193369950450545,
721        0.0194968900511231,
722        0.0196482898041878,
723        0.0197911283358190,
724        0.0199253434079123,
725        0.0200508765398072,
726        0.0201676730337687,
727        0.0202756819988200,
728        0.0203748563729175,
729        0.0204651529434560,
730        0.0205465323660984,
731        0.0206189591819181,
732        0.0206824018328499,
733        0.0207368326754401,
734        0.0207822279928917,
735        0.0208185680053983,
736        0.0208458368787627,
737        0.0208640227312962,
738        0.0208731176389954,
739        0.0208731176389954,
740        0.0208640227312962,
741        0.0208458368787627,
742        0.0208185680053983,
743        0.0207822279928917,
744        0.0207368326754401,
745        0.0206824018328499,
746        0.0206189591819181,
747        0.0205465323660984,
748        0.0204651529434560,
749        0.0203748563729175,
750        0.0202756819988200,
751        0.0201676730337687,
752        0.0200508765398072,
753        0.0199253434079123,
754        0.0197911283358190,
755        0.0196482898041878,
756        0.0194968900511231,
757        0.0193369950450545,
758        0.0191686744559934,
759        0.0189920016251754,
760        0.0188070535331042,
761        0.0186139107660094,
762        0.0184126574807331,
763        0.0182033813680609,
764        0.0179861736145128,
765        0.0177611288626114,
766        0.0175283451696437,
767        0.0172879239649355,
768        0.0170399700056559,
769        0.0167845913311726,
770        0.0165218992159766,
771        0.0162520081211971,
772        0.0159750356447283,
773        0.0156911024699895,
774        0.0154003323133401,
775        0.0151028518701744,
776        0.0147987907597169,
777        0.0144882814685445,
778        0.0141714592928592,
779        0.0138484622795371,
780        0.0135194311690004,
781        0.0131845093222756,
782        0.0128438426753249,
783        0.0124975796646449,
784        0.0121458711652067,
785        0.0117888704247183,
786        0.0114267329968529,
787        0.0110596166734735,
788        0.0106876814158841,
789        0.0103110892851360,
790        0.0099300043714212,
791        0.0095445927225849,
792        0.0091550222717888,
793        0.0087614627643580,
794        0.0083640856838475,
795        0.0079630641773633,
796        0.0075585729801782,
797        0.0071507883396855,
798        0.0067398879387430,
799        0.0063260508184704,
800        0.0059094573005900,
801        0.0054902889094487,
802        0.0050687282939456,
803        0.0046449591497966,
804        0.0042191661429919,
805        0.0037915348363451,
806        0.0033622516236779,
807        0.0029315036836558,
808        0.0024994789888943,
809        0.0020663664924131,
810        0.0016323569986067,
811        0.0011976474864367,
812        0.0007624720924706,
813        0.0003276086705538
814    ])
815)
Note: See TracBrowser for help on using the repository browser.