source: sasmodels/sasmodels/special.py @ 4073633

core_shell_microgelsmagnetic_modelticket-1257-vesicle-productticket_1156ticket_1265_superballticket_822_more_unit_tests
Last change on this file since 4073633 was 706f466, checked in by Paul Kienzle <pkienzle@…>, 7 years ago

fix doc strings for python-based sas special functions

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