source: sasmodels/sasmodels/special.py @ 110f69c

core_shell_microgelsmagnetic_modelticket-1257-vesicle-productticket_1156ticket_1265_superballticket_822_more_unit_tests
Last change on this file since 110f69c was 110f69c, checked in by Paul Kienzle <pkienzle@…>, 6 years ago

lint

  • Property mode set to 100644
File size: 20.8 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        Evaluate 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}
159            \left(1 - \frac{2!}{x^2} + \frac{4!}{x^4} - \frac{6!}{x^6} \right)
160             - \frac{\sin(x)}{x}
161            \left(\frac{1}{x} - \frac{3!}{x^3} + \frac{5!}{x^5} - \frac{7!}{x^7}\right)
162
163        For small arguments,
164
165        .. math::
166
167           \text{Si}(x) \sim x
168           - \frac{x^3}{3\times 3!} + \frac{x^5}{5 \times 5!} - \frac{x^7}{7 \times 7!}
169           + \frac{x^9}{9\times 9!} - \frac{x^{11}}{11\times 11!}
170
171    sas_3j1x_x(x):
172        Spherical Bessel form
173        $\text{sph_j1c}(x) = 3 j_1(x)/x = 3 (\sin(x) - x \cos(x))/x^3$,
174        with a limiting value of 1 at $x=0$, where $j_1(x)$ is the spherical
175        Bessel function of the first kind and first order.
176
177        This function uses a Taylor series for small $x$ for numerical accuracy.
178
179
180    sas_2J1x_x(x):
181        Bessel form $\text{sas_J1c}(x) = 2 J_1(x)/x$, with a limiting value
182        of 1 at $x=0$, where $J_1(x)$ is the Bessel function of first kind
183        and first order.
184
185
186    Gauss76Z[i], Gauss76Wt[i]:
187        Points $z_i$ and weights $w_i$ for 76-point Gaussian quadrature, respectively,
188        computing $\int_{-1}^1 f(z)\,dz \approx \sum_{i=1}^{76} w_i\,f(z_i)$.
189
190        Similar arrays are available in :code:`gauss20.c` for 20-point
191        quadrature and in :code:`gauss150.c` for 150-point quadrature.
192
193"""
194import numpy as np
195import scipy.special
196
197# Functions to add to our standard set
198from numpy import degrees, radians
199
200# C99 standard math library functions
201M_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
202from numpy import exp, log, power as pow, expm1, sqrt
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 fmin, fmax, trunc, rint
207from numpy import NAN, inf as INFINITY
208
209# erf, erfc, tgamma, lgamma  **do not use**
210
211# non-standard constants and functions
212M_PI_180, M_4PI_3 = M_PI/180, 4*M_PI/3
213
214# can't do SINCOS in python; use "s, c = SINCOS(x)" instead
215def SINCOS(x): return sin(x), cos(x)
216
217def square(x): return x*x
218
219def cube(x): return x*x*x
220
221from numpy import sinc as _sinc
222def sas_sinx_x(x): return _sinc(x/M_PI)
223def powr(x, y): return x**y
224def pown(x, n): return x**n
225
226FLOAT_SIZE = 8
227
228def polevl(x, c, n): return np.polyval(c[:n], x)
229def p1evl(x, c, n): return np.polyval(np.hstack(([1.], c))[:n], x)
230
231from scipy.special import gamma as sas_gamma
232from scipy.special import erf as sas_erf
233from scipy.special import erfc as sas_erfc
234from scipy.special import j0 as sas_J0
235from scipy.special import j1 as sas_J1
236from scipy.special import jn as sas_JN
237
238def sas_Si(x):
239    return scipy.special.sici(x)[0]
240
241def sas_j1(x):
242    if np.isscalar(x):
243        retvalue = (sin(x) - x*cos(x))/x**2 if x != 0. else 0.
244    else:
245        with np.errstate(all='ignore'):
246            retvalue = (sin(x) - x*cos(x))/x**2
247        retvalue[x == 0.] = 0.
248
249def sas_3j1x_x(x):
250    if np.isscalar(x):
251        retvalue = 3*(sin(x) - x*cos(x))/x**3 if x != 0. else 1.
252    else:
253        with np.errstate(all='ignore'):
254            retvalue = 3*(sin(x) - x*cos(x))/x**3
255        retvalue[x == 0.] = 1.
256    return retvalue
257
258def sas_2J1x_x(x):
259    if np.isscalar(x):
260        retvalue = 2*sas_J1(x)/x if x != 0 else 1.
261    else:
262        with np.errstate(all='ignore'):
263            retvalue = 2*sas_J1(x)/x
264        retvalue[x == 0] = 1.
265    return retvalue
266
267
268# Gaussians
269
270Gauss20Wt = np.array([
271    .0176140071391521,
272    .0406014298003869,
273    .0626720483341091,
274    .0832767415767047,
275    .10193011981724,
276    .118194531961518,
277    .131688638449177,
278    .142096109318382,
279    .149172986472604,
280    .152753387130726,
281    .152753387130726,
282    .149172986472604,
283    .142096109318382,
284    .131688638449177,
285    .118194531961518,
286    .10193011981724,
287    .0832767415767047,
288    .0626720483341091,
289    .0406014298003869,
290    .0176140071391521
291])
292
293Gauss20Z = np.array([
294    -.993128599185095,
295    -.963971927277914,
296    -.912234428251326,
297    -.839116971822219,
298    -.746331906460151,
299    -.636053680726515,
300    -.510867001950827,
301    -.37370608871542,
302    -.227785851141645,
303    -.076526521133497,
304    .0765265211334973,
305    .227785851141645,
306    .37370608871542,
307    .510867001950827,
308    .636053680726515,
309    .746331906460151,
310    .839116971822219,
311    .912234428251326,
312    .963971927277914,
313    .993128599185095
314])
315
316Gauss76Wt = np.array([
317    .00126779163408536,         #0
318    .00294910295364247,
319    .00462793522803742,
320    .00629918049732845,
321    .00795984747723973,
322    .00960710541471375,
323    .0112381685696677,
324    .0128502838475101,
325    .0144407317482767,
326    .0160068299122486,
327    .0175459372914742,          #10
328    .0190554584671906,
329    .020532847967908,
330    .0219756145344162,
331    .0233813253070112,
332    .0247476099206597,
333    .026072164497986,
334    .0273527555318275,
335    .028587223650054,
336    .029773487255905,
337    .0309095460374916,          #20
338    .0319934843404216,
339    .0330234743977917,
340    .0339977794120564,
341    .0349147564835508,
342    .0357728593807139,
343    .0365706411473296,
344    .0373067565423816,
345    .0379799643084053,
346    .0385891292645067,
347    .0391332242205184,          #30
348    .0396113317090621,
349    .0400226455325968,
350    .040366472122844,
351    .0406422317102947,
352    .0408494593018285,
353    .040987805464794,
354    .0410570369162294,
355    .0410570369162294,
356    .040987805464794,
357    .0408494593018285,          #40
358    .0406422317102947,
359    .040366472122844,
360    .0400226455325968,
361    .0396113317090621,
362    .0391332242205184,
363    .0385891292645067,
364    .0379799643084053,
365    .0373067565423816,
366    .0365706411473296,
367    .0357728593807139,          #50
368    .0349147564835508,
369    .0339977794120564,
370    .0330234743977917,
371    .0319934843404216,
372    .0309095460374916,
373    .029773487255905,
374    .028587223650054,
375    .0273527555318275,
376    .026072164497986,
377    .0247476099206597,          #60
378    .0233813253070112,
379    .0219756145344162,
380    .020532847967908,
381    .0190554584671906,
382    .0175459372914742,
383    .0160068299122486,
384    .0144407317482767,
385    .0128502838475101,
386    .0112381685696677,
387    .00960710541471375,         #70
388    .00795984747723973,
389    .00629918049732845,
390    .00462793522803742,
391    .00294910295364247,
392    .00126779163408536          #75 (indexed from 0)
393])
394
395Gauss76Z = np.array([
396    -.999505948362153,          #0
397    -.997397786355355,
398    -.993608772723527,
399    -.988144453359837,
400    -.981013938975656,
401    -.972229228520377,
402    -.961805126758768,
403    -.949759207710896,
404    -.936111781934811,
405    -.92088586125215,
406    -.904107119545567,          #10
407    -.885803849292083,
408    -.866006913771982,
409    -.844749694983342,
410    -.822068037328975,
411    -.7980001871612,
412    -.77258672828181,
413    -.74587051350361,
414    -.717896592387704,
415    -.688712135277641,
416    -.658366353758143,          #20
417    -.626910417672267,
418    -.594397368836793,
419    -.560882031601237,
420    -.526420920401243,
421    -.491072144462194,
422    -.454895309813726,
423    -.417951418780327,
424    -.380302767117504,
425    -.342012838966962,
426    -.303146199807908,          #30
427    -.263768387584994,
428    -.223945802196474,
429    -.183745593528914,
430    -.143235548227268,
431    -.102483975391227,
432    -.0615595913906112,
433    -.0205314039939986,
434    .0205314039939986,
435    .0615595913906112,
436    .102483975391227,                   #40
437    .143235548227268,
438    .183745593528914,
439    .223945802196474,
440    .263768387584994,
441    .303146199807908,
442    .342012838966962,
443    .380302767117504,
444    .417951418780327,
445    .454895309813726,
446    .491072144462194,           #50
447    .526420920401243,
448    .560882031601237,
449    .594397368836793,
450    .626910417672267,
451    .658366353758143,
452    .688712135277641,
453    .717896592387704,
454    .74587051350361,
455    .77258672828181,
456    .7980001871612,     #60
457    .822068037328975,
458    .844749694983342,
459    .866006913771982,
460    .885803849292083,
461    .904107119545567,
462    .92088586125215,
463    .936111781934811,
464    .949759207710896,
465    .961805126758768,
466    .972229228520377,           #70
467    .981013938975656,
468    .988144453359837,
469    .993608772723527,
470    .997397786355355,
471    .999505948362153            #75
472])
473
474Gauss150Z = np.array([
475    -0.9998723404457334,
476    -0.9993274305065947,
477    -0.9983473449340834,
478    -0.9969322929775997,
479    -0.9950828645255290,
480    -0.9927998590434373,
481    -0.9900842691660192,
482    -0.9869372772712794,
483    -0.9833602541697529,
484    -0.9793547582425894,
485    -0.9749225346595943,
486    -0.9700655145738374,
487    -0.9647858142586956,
488    -0.9590857341746905,
489    -0.9529677579610971,
490    -0.9464345513503147,
491    -0.9394889610042837,
492    -0.9321340132728527,
493    -0.9243729128743136,
494    -0.9162090414984952,
495    -0.9076459563329236,
496    -0.8986873885126239,
497    -0.8893372414942055,
498    -0.8795995893549102,
499    -0.8694786750173527,
500    -0.8589789084007133,
501    -0.8481048644991847,
502    -0.8368612813885015,
503    -0.8252530581614230,
504    -0.8132852527930605,
505    -0.8009630799369827,
506    -0.7882919086530552,
507    -0.7752772600680049,
508    -0.7619248049697269,
509    -0.7482403613363824,
510    -0.7342298918013638,
511    -0.7198995010552305,
512    -0.7052554331857488,
513    -0.6903040689571928,
514    -0.6750519230300931,
515    -0.6595056411226444,
516    -0.6436719971150083,
517    -0.6275578900977726,
518    -0.6111703413658551,
519    -0.5945164913591590,
520    -0.5776035965513142,
521    -0.5604390262878617,
522    -0.5430302595752546,
523    -0.5253848818220803,
524    -0.5075105815339176,
525    -0.4894151469632753,
526    -0.4711064627160663,
527    -0.4525925063160997,
528    -0.4338813447290861,
529    -0.4149811308476706,
530    -0.3959000999390257,
531    -0.3766465660565522,
532    -0.3572289184172501,
533    -0.3376556177463400,
534    -0.3179351925907259,
535    -0.2980762356029071,
536    -0.2780873997969574,
537    -0.2579773947782034,
538    -0.2377549829482451,
539    -0.2174289756869712,
540    -0.1970082295132342,
541    -0.1765016422258567,
542    -0.1559181490266516,
543    -0.1352667186271445,
544    -0.1145563493406956,
545    -0.0937960651617229,
546    -0.0729949118337358,
547    -0.0521619529078925,
548    -0.0313062657937972,
549    -0.0104369378042598,
550    0.0104369378042598,
551    0.0313062657937972,
552    0.0521619529078925,
553    0.0729949118337358,
554    0.0937960651617229,
555    0.1145563493406956,
556    0.1352667186271445,
557    0.1559181490266516,
558    0.1765016422258567,
559    0.1970082295132342,
560    0.2174289756869712,
561    0.2377549829482451,
562    0.2579773947782034,
563    0.2780873997969574,
564    0.2980762356029071,
565    0.3179351925907259,
566    0.3376556177463400,
567    0.3572289184172501,
568    0.3766465660565522,
569    0.3959000999390257,
570    0.4149811308476706,
571    0.4338813447290861,
572    0.4525925063160997,
573    0.4711064627160663,
574    0.4894151469632753,
575    0.5075105815339176,
576    0.5253848818220803,
577    0.5430302595752546,
578    0.5604390262878617,
579    0.5776035965513142,
580    0.5945164913591590,
581    0.6111703413658551,
582    0.6275578900977726,
583    0.6436719971150083,
584    0.6595056411226444,
585    0.6750519230300931,
586    0.6903040689571928,
587    0.7052554331857488,
588    0.7198995010552305,
589    0.7342298918013638,
590    0.7482403613363824,
591    0.7619248049697269,
592    0.7752772600680049,
593    0.7882919086530552,
594    0.8009630799369827,
595    0.8132852527930605,
596    0.8252530581614230,
597    0.8368612813885015,
598    0.8481048644991847,
599    0.8589789084007133,
600    0.8694786750173527,
601    0.8795995893549102,
602    0.8893372414942055,
603    0.8986873885126239,
604    0.9076459563329236,
605    0.9162090414984952,
606    0.9243729128743136,
607    0.9321340132728527,
608    0.9394889610042837,
609    0.9464345513503147,
610    0.9529677579610971,
611    0.9590857341746905,
612    0.9647858142586956,
613    0.9700655145738374,
614    0.9749225346595943,
615    0.9793547582425894,
616    0.9833602541697529,
617    0.9869372772712794,
618    0.9900842691660192,
619    0.9927998590434373,
620    0.9950828645255290,
621    0.9969322929775997,
622    0.9983473449340834,
623    0.9993274305065947,
624    0.9998723404457334
625])
626
627Gauss150Wt = np.array([
628    0.0003276086705538,
629    0.0007624720924706,
630    0.0011976474864367,
631    0.0016323569986067,
632    0.0020663664924131,
633    0.0024994789888943,
634    0.0029315036836558,
635    0.0033622516236779,
636    0.0037915348363451,
637    0.0042191661429919,
638    0.0046449591497966,
639    0.0050687282939456,
640    0.0054902889094487,
641    0.0059094573005900,
642    0.0063260508184704,
643    0.0067398879387430,
644    0.0071507883396855,
645    0.0075585729801782,
646    0.0079630641773633,
647    0.0083640856838475,
648    0.0087614627643580,
649    0.0091550222717888,
650    0.0095445927225849,
651    0.0099300043714212,
652    0.0103110892851360,
653    0.0106876814158841,
654    0.0110596166734735,
655    0.0114267329968529,
656    0.0117888704247183,
657    0.0121458711652067,
658    0.0124975796646449,
659    0.0128438426753249,
660    0.0131845093222756,
661    0.0135194311690004,
662    0.0138484622795371,
663    0.0141714592928592,
664    0.0144882814685445,
665    0.0147987907597169,
666    0.0151028518701744,
667    0.0154003323133401,
668    0.0156911024699895,
669    0.0159750356447283,
670    0.0162520081211971,
671    0.0165218992159766,
672    0.0167845913311726,
673    0.0170399700056559,
674    0.0172879239649355,
675    0.0175283451696437,
676    0.0177611288626114,
677    0.0179861736145128,
678    0.0182033813680609,
679    0.0184126574807331,
680    0.0186139107660094,
681    0.0188070535331042,
682    0.0189920016251754,
683    0.0191686744559934,
684    0.0193369950450545,
685    0.0194968900511231,
686    0.0196482898041878,
687    0.0197911283358190,
688    0.0199253434079123,
689    0.0200508765398072,
690    0.0201676730337687,
691    0.0202756819988200,
692    0.0203748563729175,
693    0.0204651529434560,
694    0.0205465323660984,
695    0.0206189591819181,
696    0.0206824018328499,
697    0.0207368326754401,
698    0.0207822279928917,
699    0.0208185680053983,
700    0.0208458368787627,
701    0.0208640227312962,
702    0.0208731176389954,
703    0.0208731176389954,
704    0.0208640227312962,
705    0.0208458368787627,
706    0.0208185680053983,
707    0.0207822279928917,
708    0.0207368326754401,
709    0.0206824018328499,
710    0.0206189591819181,
711    0.0205465323660984,
712    0.0204651529434560,
713    0.0203748563729175,
714    0.0202756819988200,
715    0.0201676730337687,
716    0.0200508765398072,
717    0.0199253434079123,
718    0.0197911283358190,
719    0.0196482898041878,
720    0.0194968900511231,
721    0.0193369950450545,
722    0.0191686744559934,
723    0.0189920016251754,
724    0.0188070535331042,
725    0.0186139107660094,
726    0.0184126574807331,
727    0.0182033813680609,
728    0.0179861736145128,
729    0.0177611288626114,
730    0.0175283451696437,
731    0.0172879239649355,
732    0.0170399700056559,
733    0.0167845913311726,
734    0.0165218992159766,
735    0.0162520081211971,
736    0.0159750356447283,
737    0.0156911024699895,
738    0.0154003323133401,
739    0.0151028518701744,
740    0.0147987907597169,
741    0.0144882814685445,
742    0.0141714592928592,
743    0.0138484622795371,
744    0.0135194311690004,
745    0.0131845093222756,
746    0.0128438426753249,
747    0.0124975796646449,
748    0.0121458711652067,
749    0.0117888704247183,
750    0.0114267329968529,
751    0.0110596166734735,
752    0.0106876814158841,
753    0.0103110892851360,
754    0.0099300043714212,
755    0.0095445927225849,
756    0.0091550222717888,
757    0.0087614627643580,
758    0.0083640856838475,
759    0.0079630641773633,
760    0.0075585729801782,
761    0.0071507883396855,
762    0.0067398879387430,
763    0.0063260508184704,
764    0.0059094573005900,
765    0.0054902889094487,
766    0.0050687282939456,
767    0.0046449591497966,
768    0.0042191661429919,
769    0.0037915348363451,
770    0.0033622516236779,
771    0.0029315036836558,
772    0.0024994789888943,
773    0.0020663664924131,
774    0.0016323569986067,
775    0.0011976474864367,
776    0.0007624720924706,
777    0.0003276086705538
778])
Note: See TracBrowser for help on using the repository browser.