1 /* Math module -- standard C math library functions, pi and e */
2
3 /* Here are some comments from Tim Peters, extracted from the
4 discussion attached to http://bugs.python.org/issue1640. They
5 describe the general aims of the math module with respect to
6 special values, IEEE-754 floating-point exceptions, and Python
7 exceptions.
8
9 These are the "spirit of 754" rules:
10
11 1. If the mathematical result is a real number, but of magnitude too
12 large to approximate by a machine float, overflow is signaled and the
13 result is an infinity (with the appropriate sign).
14
15 2. If the mathematical result is a real number, but of magnitude too
16 small to approximate by a machine float, underflow is signaled and the
17 result is a zero (with the appropriate sign).
18
19 3. At a singularity (a value x such that the limit of f(y) as y
20 approaches x exists and is an infinity), "divide by zero" is signaled
21 and the result is an infinity (with the appropriate sign). This is
22 complicated a little by that the left-side and right-side limits may
23 not be the same; e.g., 1/x approaches +inf or -inf as x approaches 0
24 from the positive or negative directions. In that specific case, the
25 sign of the zero determines the result of 1/0.
26
27 4. At a point where a function has no defined result in the extended
28 reals (i.e., the reals plus an infinity or two), invalid operation is
29 signaled and a NaN is returned.
30
31 And these are what Python has historically /tried/ to do (but not
32 always successfully, as platform libm behavior varies a lot):
33
34 For #1, raise OverflowError.
35
36 For #2, return a zero (with the appropriate sign if that happens by
37 accident ;-)).
38
39 For #3 and #4, raise ValueError. It may have made sense to raise
40 Python's ZeroDivisionError in #3, but historically that's only been
41 raised for division by zero and mod by zero.
42
43 */
44
45 /*
46 In general, on an IEEE-754 platform the aim is to follow the C99
47 standard, including Annex 'F', whenever possible. Where the
48 standard recommends raising the 'divide-by-zero' or 'invalid'
49 floating-point exceptions, Python should raise a ValueError. Where
50 the standard recommends raising 'overflow', Python should raise an
51 OverflowError. In all other circumstances a value should be
52 returned.
53 */
54
55 #include "Python.h"
56 #include "_math.h"
57
58 #include "clinic/mathmodule.c.h"
59
60 /*[clinic input]
61 module math
62 [clinic start generated code]*/
63 /*[clinic end generated code: output=da39a3ee5e6b4b0d input=76bc7002685dd942]*/
64
65
66 /*
67 sin(pi*x), giving accurate results for all finite x (especially x
68 integral or close to an integer). This is here for use in the
69 reflection formula for the gamma function. It conforms to IEEE
70 754-2008 for finite arguments, but not for infinities or nans.
71 */
72
73 static const double pi = 3.141592653589793238462643383279502884197;
74 static const double logpi = 1.144729885849400174143427351353058711647;
75 #if !defined(HAVE_ERF) || !defined(HAVE_ERFC)
76 static const double sqrtpi = 1.772453850905516027298167483341145182798;
77 #endif /* !defined(HAVE_ERF) || !defined(HAVE_ERFC) */
78
79
80 /* Version of PyFloat_AsDouble() with in-line fast paths
81 for exact floats and integers. Gives a substantial
82 speed improvement for extracting float arguments.
83 */
84
85 #define ASSIGN_DOUBLE(target_var, obj, error_label) \
86 if (PyFloat_CheckExact(obj)) { \
87 target_var = PyFloat_AS_DOUBLE(obj); \
88 } \
89 else if (PyLong_CheckExact(obj)) { \
90 target_var = PyLong_AsDouble(obj); \
91 if (target_var == -1.0 && PyErr_Occurred()) { \
92 goto error_label; \
93 } \
94 } \
95 else { \
96 target_var = PyFloat_AsDouble(obj); \
97 if (target_var == -1.0 && PyErr_Occurred()) { \
98 goto error_label; \
99 } \
100 }
101
102 static double
m_sinpi(double x)103 m_sinpi(double x)
104 {
105 double y, r;
106 int n;
107 /* this function should only ever be called for finite arguments */
108 assert(Py_IS_FINITE(x));
109 y = fmod(fabs(x), 2.0);
110 n = (int)round(2.0*y);
111 assert(0 <= n && n <= 4);
112 switch (n) {
113 case 0:
114 r = sin(pi*y);
115 break;
116 case 1:
117 r = cos(pi*(y-0.5));
118 break;
119 case 2:
120 /* N.B. -sin(pi*(y-1.0)) is *not* equivalent: it would give
121 -0.0 instead of 0.0 when y == 1.0. */
122 r = sin(pi*(1.0-y));
123 break;
124 case 3:
125 r = -cos(pi*(y-1.5));
126 break;
127 case 4:
128 r = sin(pi*(y-2.0));
129 break;
130 default:
131 Py_UNREACHABLE();
132 }
133 return copysign(1.0, x)*r;
134 }
135
136 /* Implementation of the real gamma function. In extensive but non-exhaustive
137 random tests, this function proved accurate to within <= 10 ulps across the
138 entire float domain. Note that accuracy may depend on the quality of the
139 system math functions, the pow function in particular. Special cases
140 follow C99 annex F. The parameters and method are tailored to platforms
141 whose double format is the IEEE 754 binary64 format.
142
143 Method: for x > 0.0 we use the Lanczos approximation with parameters N=13
144 and g=6.024680040776729583740234375; these parameters are amongst those
145 used by the Boost library. Following Boost (again), we re-express the
146 Lanczos sum as a rational function, and compute it that way. The
147 coefficients below were computed independently using MPFR, and have been
148 double-checked against the coefficients in the Boost source code.
149
150 For x < 0.0 we use the reflection formula.
151
152 There's one minor tweak that deserves explanation: Lanczos' formula for
153 Gamma(x) involves computing pow(x+g-0.5, x-0.5) / exp(x+g-0.5). For many x
154 values, x+g-0.5 can be represented exactly. However, in cases where it
155 can't be represented exactly the small error in x+g-0.5 can be magnified
156 significantly by the pow and exp calls, especially for large x. A cheap
157 correction is to multiply by (1 + e*g/(x+g-0.5)), where e is the error
158 involved in the computation of x+g-0.5 (that is, e = computed value of
159 x+g-0.5 - exact value of x+g-0.5). Here's the proof:
160
161 Correction factor
162 -----------------
163 Write x+g-0.5 = y-e, where y is exactly representable as an IEEE 754
164 double, and e is tiny. Then:
165
166 pow(x+g-0.5,x-0.5)/exp(x+g-0.5) = pow(y-e, x-0.5)/exp(y-e)
167 = pow(y, x-0.5)/exp(y) * C,
168
169 where the correction_factor C is given by
170
171 C = pow(1-e/y, x-0.5) * exp(e)
172
173 Since e is tiny, pow(1-e/y, x-0.5) ~ 1-(x-0.5)*e/y, and exp(x) ~ 1+e, so:
174
175 C ~ (1-(x-0.5)*e/y) * (1+e) ~ 1 + e*(y-(x-0.5))/y
176
177 But y-(x-0.5) = g+e, and g+e ~ g. So we get C ~ 1 + e*g/y, and
178
179 pow(x+g-0.5,x-0.5)/exp(x+g-0.5) ~ pow(y, x-0.5)/exp(y) * (1 + e*g/y),
180
181 Note that for accuracy, when computing r*C it's better to do
182
183 r + e*g/y*r;
184
185 than
186
187 r * (1 + e*g/y);
188
189 since the addition in the latter throws away most of the bits of
190 information in e*g/y.
191 */
192
193 #define LANCZOS_N 13
194 static const double lanczos_g = 6.024680040776729583740234375;
195 static const double lanczos_g_minus_half = 5.524680040776729583740234375;
196 static const double lanczos_num_coeffs[LANCZOS_N] = {
197 23531376880.410759688572007674451636754734846804940,
198 42919803642.649098768957899047001988850926355848959,
199 35711959237.355668049440185451547166705960488635843,
200 17921034426.037209699919755754458931112671403265390,
201 6039542586.3520280050642916443072979210699388420708,
202 1439720407.3117216736632230727949123939715485786772,
203 248874557.86205415651146038641322942321632125127801,
204 31426415.585400194380614231628318205362874684987640,
205 2876370.6289353724412254090516208496135991145378768,
206 186056.26539522349504029498971604569928220784236328,
207 8071.6720023658162106380029022722506138218516325024,
208 210.82427775157934587250973392071336271166969580291,
209 2.5066282746310002701649081771338373386264310793408
210 };
211
212 /* denominator is x*(x+1)*...*(x+LANCZOS_N-2) */
213 static const double lanczos_den_coeffs[LANCZOS_N] = {
214 0.0, 39916800.0, 120543840.0, 150917976.0, 105258076.0, 45995730.0,
215 13339535.0, 2637558.0, 357423.0, 32670.0, 1925.0, 66.0, 1.0};
216
217 /* gamma values for small positive integers, 1 though NGAMMA_INTEGRAL */
218 #define NGAMMA_INTEGRAL 23
219 static const double gamma_integral[NGAMMA_INTEGRAL] = {
220 1.0, 1.0, 2.0, 6.0, 24.0, 120.0, 720.0, 5040.0, 40320.0, 362880.0,
221 3628800.0, 39916800.0, 479001600.0, 6227020800.0, 87178291200.0,
222 1307674368000.0, 20922789888000.0, 355687428096000.0,
223 6402373705728000.0, 121645100408832000.0, 2432902008176640000.0,
224 51090942171709440000.0, 1124000727777607680000.0,
225 };
226
227 /* Lanczos' sum L_g(x), for positive x */
228
229 static double
lanczos_sum(double x)230 lanczos_sum(double x)
231 {
232 double num = 0.0, den = 0.0;
233 int i;
234 assert(x > 0.0);
235 /* evaluate the rational function lanczos_sum(x). For large
236 x, the obvious algorithm risks overflow, so we instead
237 rescale the denominator and numerator of the rational
238 function by x**(1-LANCZOS_N) and treat this as a
239 rational function in 1/x. This also reduces the error for
240 larger x values. The choice of cutoff point (5.0 below) is
241 somewhat arbitrary; in tests, smaller cutoff values than
242 this resulted in lower accuracy. */
243 if (x < 5.0) {
244 for (i = LANCZOS_N; --i >= 0; ) {
245 num = num * x + lanczos_num_coeffs[i];
246 den = den * x + lanczos_den_coeffs[i];
247 }
248 }
249 else {
250 for (i = 0; i < LANCZOS_N; i++) {
251 num = num / x + lanczos_num_coeffs[i];
252 den = den / x + lanczos_den_coeffs[i];
253 }
254 }
255 return num/den;
256 }
257
258 /* Constant for +infinity, generated in the same way as float('inf'). */
259
260 static double
m_inf(void)261 m_inf(void)
262 {
263 #ifndef PY_NO_SHORT_FLOAT_REPR
264 return _Py_dg_infinity(0);
265 #else
266 return Py_HUGE_VAL;
267 #endif
268 }
269
270 /* Constant nan value, generated in the same way as float('nan'). */
271 /* We don't currently assume that Py_NAN is defined everywhere. */
272
273 #if !defined(PY_NO_SHORT_FLOAT_REPR) || defined(Py_NAN)
274
275 static double
m_nan(void)276 m_nan(void)
277 {
278 #ifndef PY_NO_SHORT_FLOAT_REPR
279 return _Py_dg_stdnan(0);
280 #else
281 return Py_NAN;
282 #endif
283 }
284
285 #endif
286
287 static double
m_tgamma(double x)288 m_tgamma(double x)
289 {
290 double absx, r, y, z, sqrtpow;
291
292 /* special cases */
293 if (!Py_IS_FINITE(x)) {
294 if (Py_IS_NAN(x) || x > 0.0)
295 return x; /* tgamma(nan) = nan, tgamma(inf) = inf */
296 else {
297 errno = EDOM;
298 return Py_NAN; /* tgamma(-inf) = nan, invalid */
299 }
300 }
301 if (x == 0.0) {
302 errno = EDOM;
303 /* tgamma(+-0.0) = +-inf, divide-by-zero */
304 return copysign(Py_HUGE_VAL, x);
305 }
306
307 /* integer arguments */
308 if (x == floor(x)) {
309 if (x < 0.0) {
310 errno = EDOM; /* tgamma(n) = nan, invalid for */
311 return Py_NAN; /* negative integers n */
312 }
313 if (x <= NGAMMA_INTEGRAL)
314 return gamma_integral[(int)x - 1];
315 }
316 absx = fabs(x);
317
318 /* tiny arguments: tgamma(x) ~ 1/x for x near 0 */
319 if (absx < 1e-20) {
320 r = 1.0/x;
321 if (Py_IS_INFINITY(r))
322 errno = ERANGE;
323 return r;
324 }
325
326 /* large arguments: assuming IEEE 754 doubles, tgamma(x) overflows for
327 x > 200, and underflows to +-0.0 for x < -200, not a negative
328 integer. */
329 if (absx > 200.0) {
330 if (x < 0.0) {
331 return 0.0/m_sinpi(x);
332 }
333 else {
334 errno = ERANGE;
335 return Py_HUGE_VAL;
336 }
337 }
338
339 y = absx + lanczos_g_minus_half;
340 /* compute error in sum */
341 if (absx > lanczos_g_minus_half) {
342 /* note: the correction can be foiled by an optimizing
343 compiler that (incorrectly) thinks that an expression like
344 a + b - a - b can be optimized to 0.0. This shouldn't
345 happen in a standards-conforming compiler. */
346 double q = y - absx;
347 z = q - lanczos_g_minus_half;
348 }
349 else {
350 double q = y - lanczos_g_minus_half;
351 z = q - absx;
352 }
353 z = z * lanczos_g / y;
354 if (x < 0.0) {
355 r = -pi / m_sinpi(absx) / absx * exp(y) / lanczos_sum(absx);
356 r -= z * r;
357 if (absx < 140.0) {
358 r /= pow(y, absx - 0.5);
359 }
360 else {
361 sqrtpow = pow(y, absx / 2.0 - 0.25);
362 r /= sqrtpow;
363 r /= sqrtpow;
364 }
365 }
366 else {
367 r = lanczos_sum(absx) / exp(y);
368 r += z * r;
369 if (absx < 140.0) {
370 r *= pow(y, absx - 0.5);
371 }
372 else {
373 sqrtpow = pow(y, absx / 2.0 - 0.25);
374 r *= sqrtpow;
375 r *= sqrtpow;
376 }
377 }
378 if (Py_IS_INFINITY(r))
379 errno = ERANGE;
380 return r;
381 }
382
383 /*
384 lgamma: natural log of the absolute value of the Gamma function.
385 For large arguments, Lanczos' formula works extremely well here.
386 */
387
388 static double
m_lgamma(double x)389 m_lgamma(double x)
390 {
391 double r;
392 double absx;
393
394 /* special cases */
395 if (!Py_IS_FINITE(x)) {
396 if (Py_IS_NAN(x))
397 return x; /* lgamma(nan) = nan */
398 else
399 return Py_HUGE_VAL; /* lgamma(+-inf) = +inf */
400 }
401
402 /* integer arguments */
403 if (x == floor(x) && x <= 2.0) {
404 if (x <= 0.0) {
405 errno = EDOM; /* lgamma(n) = inf, divide-by-zero for */
406 return Py_HUGE_VAL; /* integers n <= 0 */
407 }
408 else {
409 return 0.0; /* lgamma(1) = lgamma(2) = 0.0 */
410 }
411 }
412
413 absx = fabs(x);
414 /* tiny arguments: lgamma(x) ~ -log(fabs(x)) for small x */
415 if (absx < 1e-20)
416 return -log(absx);
417
418 /* Lanczos' formula. We could save a fraction of a ulp in accuracy by
419 having a second set of numerator coefficients for lanczos_sum that
420 absorbed the exp(-lanczos_g) term, and throwing out the lanczos_g
421 subtraction below; it's probably not worth it. */
422 r = log(lanczos_sum(absx)) - lanczos_g;
423 r += (absx - 0.5) * (log(absx + lanczos_g - 0.5) - 1);
424 if (x < 0.0)
425 /* Use reflection formula to get value for negative x. */
426 r = logpi - log(fabs(m_sinpi(absx))) - log(absx) - r;
427 if (Py_IS_INFINITY(r))
428 errno = ERANGE;
429 return r;
430 }
431
432 #if !defined(HAVE_ERF) || !defined(HAVE_ERFC)
433
434 /*
435 Implementations of the error function erf(x) and the complementary error
436 function erfc(x).
437
438 Method: we use a series approximation for erf for small x, and a continued
439 fraction approximation for erfc(x) for larger x;
440 combined with the relations erf(-x) = -erf(x) and erfc(x) = 1.0 - erf(x),
441 this gives us erf(x) and erfc(x) for all x.
442
443 The series expansion used is:
444
445 erf(x) = x*exp(-x*x)/sqrt(pi) * [
446 2/1 + 4/3 x**2 + 8/15 x**4 + 16/105 x**6 + ...]
447
448 The coefficient of x**(2k-2) here is 4**k*factorial(k)/factorial(2*k).
449 This series converges well for smallish x, but slowly for larger x.
450
451 The continued fraction expansion used is:
452
453 erfc(x) = x*exp(-x*x)/sqrt(pi) * [1/(0.5 + x**2 -) 0.5/(2.5 + x**2 - )
454 3.0/(4.5 + x**2 - ) 7.5/(6.5 + x**2 - ) ...]
455
456 after the first term, the general term has the form:
457
458 k*(k-0.5)/(2*k+0.5 + x**2 - ...).
459
460 This expansion converges fast for larger x, but convergence becomes
461 infinitely slow as x approaches 0.0. The (somewhat naive) continued
462 fraction evaluation algorithm used below also risks overflow for large x;
463 but for large x, erfc(x) == 0.0 to within machine precision. (For
464 example, erfc(30.0) is approximately 2.56e-393).
465
466 Parameters: use series expansion for abs(x) < ERF_SERIES_CUTOFF and
467 continued fraction expansion for ERF_SERIES_CUTOFF <= abs(x) <
468 ERFC_CONTFRAC_CUTOFF. ERFC_SERIES_TERMS and ERFC_CONTFRAC_TERMS are the
469 numbers of terms to use for the relevant expansions. */
470
471 #define ERF_SERIES_CUTOFF 1.5
472 #define ERF_SERIES_TERMS 25
473 #define ERFC_CONTFRAC_CUTOFF 30.0
474 #define ERFC_CONTFRAC_TERMS 50
475
476 /*
477 Error function, via power series.
478
479 Given a finite float x, return an approximation to erf(x).
480 Converges reasonably fast for small x.
481 */
482
483 static double
m_erf_series(double x)484 m_erf_series(double x)
485 {
486 double x2, acc, fk, result;
487 int i, saved_errno;
488
489 x2 = x * x;
490 acc = 0.0;
491 fk = (double)ERF_SERIES_TERMS + 0.5;
492 for (i = 0; i < ERF_SERIES_TERMS; i++) {
493 acc = 2.0 + x2 * acc / fk;
494 fk -= 1.0;
495 }
496 /* Make sure the exp call doesn't affect errno;
497 see m_erfc_contfrac for more. */
498 saved_errno = errno;
499 result = acc * x * exp(-x2) / sqrtpi;
500 errno = saved_errno;
501 return result;
502 }
503
504 /*
505 Complementary error function, via continued fraction expansion.
506
507 Given a positive float x, return an approximation to erfc(x). Converges
508 reasonably fast for x large (say, x > 2.0), and should be safe from
509 overflow if x and nterms are not too large. On an IEEE 754 machine, with x
510 <= 30.0, we're safe up to nterms = 100. For x >= 30.0, erfc(x) is smaller
511 than the smallest representable nonzero float. */
512
513 static double
m_erfc_contfrac(double x)514 m_erfc_contfrac(double x)
515 {
516 double x2, a, da, p, p_last, q, q_last, b, result;
517 int i, saved_errno;
518
519 if (x >= ERFC_CONTFRAC_CUTOFF)
520 return 0.0;
521
522 x2 = x*x;
523 a = 0.0;
524 da = 0.5;
525 p = 1.0; p_last = 0.0;
526 q = da + x2; q_last = 1.0;
527 for (i = 0; i < ERFC_CONTFRAC_TERMS; i++) {
528 double temp;
529 a += da;
530 da += 2.0;
531 b = da + x2;
532 temp = p; p = b*p - a*p_last; p_last = temp;
533 temp = q; q = b*q - a*q_last; q_last = temp;
534 }
535 /* Issue #8986: On some platforms, exp sets errno on underflow to zero;
536 save the current errno value so that we can restore it later. */
537 saved_errno = errno;
538 result = p / q * x * exp(-x2) / sqrtpi;
539 errno = saved_errno;
540 return result;
541 }
542
543 #endif /* !defined(HAVE_ERF) || !defined(HAVE_ERFC) */
544
545 /* Error function erf(x), for general x */
546
547 static double
m_erf(double x)548 m_erf(double x)
549 {
550 #ifdef HAVE_ERF
551 return erf(x);
552 #else
553 double absx, cf;
554
555 if (Py_IS_NAN(x))
556 return x;
557 absx = fabs(x);
558 if (absx < ERF_SERIES_CUTOFF)
559 return m_erf_series(x);
560 else {
561 cf = m_erfc_contfrac(absx);
562 return x > 0.0 ? 1.0 - cf : cf - 1.0;
563 }
564 #endif
565 }
566
567 /* Complementary error function erfc(x), for general x. */
568
569 static double
m_erfc(double x)570 m_erfc(double x)
571 {
572 #ifdef HAVE_ERFC
573 return erfc(x);
574 #else
575 double absx, cf;
576
577 if (Py_IS_NAN(x))
578 return x;
579 absx = fabs(x);
580 if (absx < ERF_SERIES_CUTOFF)
581 return 1.0 - m_erf_series(x);
582 else {
583 cf = m_erfc_contfrac(absx);
584 return x > 0.0 ? cf : 2.0 - cf;
585 }
586 #endif
587 }
588
589 /*
590 wrapper for atan2 that deals directly with special cases before
591 delegating to the platform libm for the remaining cases. This
592 is necessary to get consistent behaviour across platforms.
593 Windows, FreeBSD and alpha Tru64 are amongst platforms that don't
594 always follow C99.
595 */
596
597 static double
m_atan2(double y,double x)598 m_atan2(double y, double x)
599 {
600 if (Py_IS_NAN(x) || Py_IS_NAN(y))
601 return Py_NAN;
602 if (Py_IS_INFINITY(y)) {
603 if (Py_IS_INFINITY(x)) {
604 if (copysign(1., x) == 1.)
605 /* atan2(+-inf, +inf) == +-pi/4 */
606 return copysign(0.25*Py_MATH_PI, y);
607 else
608 /* atan2(+-inf, -inf) == +-pi*3/4 */
609 return copysign(0.75*Py_MATH_PI, y);
610 }
611 /* atan2(+-inf, x) == +-pi/2 for finite x */
612 return copysign(0.5*Py_MATH_PI, y);
613 }
614 if (Py_IS_INFINITY(x) || y == 0.) {
615 if (copysign(1., x) == 1.)
616 /* atan2(+-y, +inf) = atan2(+-0, +x) = +-0. */
617 return copysign(0., y);
618 else
619 /* atan2(+-y, -inf) = atan2(+-0., -x) = +-pi. */
620 return copysign(Py_MATH_PI, y);
621 }
622 return atan2(y, x);
623 }
624
625
626 /* IEEE 754-style remainder operation: x - n*y where n*y is the nearest
627 multiple of y to x, taking n even in the case of a tie. Assuming an IEEE 754
628 binary floating-point format, the result is always exact. */
629
630 static double
m_remainder(double x,double y)631 m_remainder(double x, double y)
632 {
633 /* Deal with most common case first. */
634 if (Py_IS_FINITE(x) && Py_IS_FINITE(y)) {
635 double absx, absy, c, m, r;
636
637 if (y == 0.0) {
638 return Py_NAN;
639 }
640
641 absx = fabs(x);
642 absy = fabs(y);
643 m = fmod(absx, absy);
644
645 /*
646 Warning: some subtlety here. What we *want* to know at this point is
647 whether the remainder m is less than, equal to, or greater than half
648 of absy. However, we can't do that comparison directly because we
649 can't be sure that 0.5*absy is representable (the mutiplication
650 might incur precision loss due to underflow). So instead we compare
651 m with the complement c = absy - m: m < 0.5*absy if and only if m <
652 c, and so on. The catch is that absy - m might also not be
653 representable, but it turns out that it doesn't matter:
654
655 - if m > 0.5*absy then absy - m is exactly representable, by
656 Sterbenz's lemma, so m > c
657 - if m == 0.5*absy then again absy - m is exactly representable
658 and m == c
659 - if m < 0.5*absy then either (i) 0.5*absy is exactly representable,
660 in which case 0.5*absy < absy - m, so 0.5*absy <= c and hence m <
661 c, or (ii) absy is tiny, either subnormal or in the lowest normal
662 binade. Then absy - m is exactly representable and again m < c.
663 */
664
665 c = absy - m;
666 if (m < c) {
667 r = m;
668 }
669 else if (m > c) {
670 r = -c;
671 }
672 else {
673 /*
674 Here absx is exactly halfway between two multiples of absy,
675 and we need to choose the even multiple. x now has the form
676
677 absx = n * absy + m
678
679 for some integer n (recalling that m = 0.5*absy at this point).
680 If n is even we want to return m; if n is odd, we need to
681 return -m.
682
683 So
684
685 0.5 * (absx - m) = (n/2) * absy
686
687 and now reducing modulo absy gives us:
688
689 | m, if n is odd
690 fmod(0.5 * (absx - m), absy) = |
691 | 0, if n is even
692
693 Now m - 2.0 * fmod(...) gives the desired result: m
694 if n is even, -m if m is odd.
695
696 Note that all steps in fmod(0.5 * (absx - m), absy)
697 will be computed exactly, with no rounding error
698 introduced.
699 */
700 assert(m == c);
701 r = m - 2.0 * fmod(0.5 * (absx - m), absy);
702 }
703 return copysign(1.0, x) * r;
704 }
705
706 /* Special values. */
707 if (Py_IS_NAN(x)) {
708 return x;
709 }
710 if (Py_IS_NAN(y)) {
711 return y;
712 }
713 if (Py_IS_INFINITY(x)) {
714 return Py_NAN;
715 }
716 assert(Py_IS_INFINITY(y));
717 return x;
718 }
719
720
721 /*
722 Various platforms (Solaris, OpenBSD) do nonstandard things for log(0),
723 log(-ve), log(NaN). Here are wrappers for log and log10 that deal with
724 special values directly, passing positive non-special values through to
725 the system log/log10.
726 */
727
728 static double
m_log(double x)729 m_log(double x)
730 {
731 if (Py_IS_FINITE(x)) {
732 if (x > 0.0)
733 return log(x);
734 errno = EDOM;
735 if (x == 0.0)
736 return -Py_HUGE_VAL; /* log(0) = -inf */
737 else
738 return Py_NAN; /* log(-ve) = nan */
739 }
740 else if (Py_IS_NAN(x))
741 return x; /* log(nan) = nan */
742 else if (x > 0.0)
743 return x; /* log(inf) = inf */
744 else {
745 errno = EDOM;
746 return Py_NAN; /* log(-inf) = nan */
747 }
748 }
749
750 /*
751 log2: log to base 2.
752
753 Uses an algorithm that should:
754
755 (a) produce exact results for powers of 2, and
756 (b) give a monotonic log2 (for positive finite floats),
757 assuming that the system log is monotonic.
758 */
759
760 static double
m_log2(double x)761 m_log2(double x)
762 {
763 if (!Py_IS_FINITE(x)) {
764 if (Py_IS_NAN(x))
765 return x; /* log2(nan) = nan */
766 else if (x > 0.0)
767 return x; /* log2(+inf) = +inf */
768 else {
769 errno = EDOM;
770 return Py_NAN; /* log2(-inf) = nan, invalid-operation */
771 }
772 }
773
774 if (x > 0.0) {
775 #ifdef HAVE_LOG2
776 return log2(x);
777 #else
778 double m;
779 int e;
780 m = frexp(x, &e);
781 /* We want log2(m * 2**e) == log(m) / log(2) + e. Care is needed when
782 * x is just greater than 1.0: in that case e is 1, log(m) is negative,
783 * and we get significant cancellation error from the addition of
784 * log(m) / log(2) to e. The slight rewrite of the expression below
785 * avoids this problem.
786 */
787 if (x >= 1.0) {
788 return log(2.0 * m) / log(2.0) + (e - 1);
789 }
790 else {
791 return log(m) / log(2.0) + e;
792 }
793 #endif
794 }
795 else if (x == 0.0) {
796 errno = EDOM;
797 return -Py_HUGE_VAL; /* log2(0) = -inf, divide-by-zero */
798 }
799 else {
800 errno = EDOM;
801 return Py_NAN; /* log2(-inf) = nan, invalid-operation */
802 }
803 }
804
805 static double
m_log10(double x)806 m_log10(double x)
807 {
808 if (Py_IS_FINITE(x)) {
809 if (x > 0.0)
810 return log10(x);
811 errno = EDOM;
812 if (x == 0.0)
813 return -Py_HUGE_VAL; /* log10(0) = -inf */
814 else
815 return Py_NAN; /* log10(-ve) = nan */
816 }
817 else if (Py_IS_NAN(x))
818 return x; /* log10(nan) = nan */
819 else if (x > 0.0)
820 return x; /* log10(inf) = inf */
821 else {
822 errno = EDOM;
823 return Py_NAN; /* log10(-inf) = nan */
824 }
825 }
826
827
828 /*[clinic input]
829 math.gcd
830
831 x as a: object
832 y as b: object
833 /
834
835 greatest common divisor of x and y
836 [clinic start generated code]*/
837
838 static PyObject *
math_gcd_impl(PyObject * module,PyObject * a,PyObject * b)839 math_gcd_impl(PyObject *module, PyObject *a, PyObject *b)
840 /*[clinic end generated code: output=7b2e0c151bd7a5d8 input=c2691e57fb2a98fa]*/
841 {
842 PyObject *g;
843
844 a = PyNumber_Index(a);
845 if (a == NULL)
846 return NULL;
847 b = PyNumber_Index(b);
848 if (b == NULL) {
849 Py_DECREF(a);
850 return NULL;
851 }
852 g = _PyLong_GCD(a, b);
853 Py_DECREF(a);
854 Py_DECREF(b);
855 return g;
856 }
857
858
859 /* Call is_error when errno != 0, and where x is the result libm
860 * returned. is_error will usually set up an exception and return
861 * true (1), but may return false (0) without setting up an exception.
862 */
863 static int
is_error(double x)864 is_error(double x)
865 {
866 int result = 1; /* presumption of guilt */
867 assert(errno); /* non-zero errno is a precondition for calling */
868 if (errno == EDOM)
869 PyErr_SetString(PyExc_ValueError, "math domain error");
870
871 else if (errno == ERANGE) {
872 /* ANSI C generally requires libm functions to set ERANGE
873 * on overflow, but also generally *allows* them to set
874 * ERANGE on underflow too. There's no consistency about
875 * the latter across platforms.
876 * Alas, C99 never requires that errno be set.
877 * Here we suppress the underflow errors (libm functions
878 * should return a zero on underflow, and +- HUGE_VAL on
879 * overflow, so testing the result for zero suffices to
880 * distinguish the cases).
881 *
882 * On some platforms (Ubuntu/ia64) it seems that errno can be
883 * set to ERANGE for subnormal results that do *not* underflow
884 * to zero. So to be safe, we'll ignore ERANGE whenever the
885 * function result is less than one in absolute value.
886 */
887 if (fabs(x) < 1.0)
888 result = 0;
889 else
890 PyErr_SetString(PyExc_OverflowError,
891 "math range error");
892 }
893 else
894 /* Unexpected math error */
895 PyErr_SetFromErrno(PyExc_ValueError);
896 return result;
897 }
898
899 /*
900 math_1 is used to wrap a libm function f that takes a double
901 argument and returns a double.
902
903 The error reporting follows these rules, which are designed to do
904 the right thing on C89/C99 platforms and IEEE 754/non IEEE 754
905 platforms.
906
907 - a NaN result from non-NaN inputs causes ValueError to be raised
908 - an infinite result from finite inputs causes OverflowError to be
909 raised if can_overflow is 1, or raises ValueError if can_overflow
910 is 0.
911 - if the result is finite and errno == EDOM then ValueError is
912 raised
913 - if the result is finite and nonzero and errno == ERANGE then
914 OverflowError is raised
915
916 The last rule is used to catch overflow on platforms which follow
917 C89 but for which HUGE_VAL is not an infinity.
918
919 For the majority of one-argument functions these rules are enough
920 to ensure that Python's functions behave as specified in 'Annex F'
921 of the C99 standard, with the 'invalid' and 'divide-by-zero'
922 floating-point exceptions mapping to Python's ValueError and the
923 'overflow' floating-point exception mapping to OverflowError.
924 math_1 only works for functions that don't have singularities *and*
925 the possibility of overflow; fortunately, that covers everything we
926 care about right now.
927 */
928
929 static PyObject *
math_1_to_whatever(PyObject * arg,double (* func)(double),PyObject * (* from_double_func)(double),int can_overflow)930 math_1_to_whatever(PyObject *arg, double (*func) (double),
931 PyObject *(*from_double_func) (double),
932 int can_overflow)
933 {
934 double x, r;
935 x = PyFloat_AsDouble(arg);
936 if (x == -1.0 && PyErr_Occurred())
937 return NULL;
938 errno = 0;
939 PyFPE_START_PROTECT("in math_1", return 0);
940 r = (*func)(x);
941 PyFPE_END_PROTECT(r);
942 if (Py_IS_NAN(r) && !Py_IS_NAN(x)) {
943 PyErr_SetString(PyExc_ValueError,
944 "math domain error"); /* invalid arg */
945 return NULL;
946 }
947 if (Py_IS_INFINITY(r) && Py_IS_FINITE(x)) {
948 if (can_overflow)
949 PyErr_SetString(PyExc_OverflowError,
950 "math range error"); /* overflow */
951 else
952 PyErr_SetString(PyExc_ValueError,
953 "math domain error"); /* singularity */
954 return NULL;
955 }
956 if (Py_IS_FINITE(r) && errno && is_error(r))
957 /* this branch unnecessary on most platforms */
958 return NULL;
959
960 return (*from_double_func)(r);
961 }
962
963 /* variant of math_1, to be used when the function being wrapped is known to
964 set errno properly (that is, errno = EDOM for invalid or divide-by-zero,
965 errno = ERANGE for overflow). */
966
967 static PyObject *
math_1a(PyObject * arg,double (* func)(double))968 math_1a(PyObject *arg, double (*func) (double))
969 {
970 double x, r;
971 x = PyFloat_AsDouble(arg);
972 if (x == -1.0 && PyErr_Occurred())
973 return NULL;
974 errno = 0;
975 PyFPE_START_PROTECT("in math_1a", return 0);
976 r = (*func)(x);
977 PyFPE_END_PROTECT(r);
978 if (errno && is_error(r))
979 return NULL;
980 return PyFloat_FromDouble(r);
981 }
982
983 /*
984 math_2 is used to wrap a libm function f that takes two double
985 arguments and returns a double.
986
987 The error reporting follows these rules, which are designed to do
988 the right thing on C89/C99 platforms and IEEE 754/non IEEE 754
989 platforms.
990
991 - a NaN result from non-NaN inputs causes ValueError to be raised
992 - an infinite result from finite inputs causes OverflowError to be
993 raised.
994 - if the result is finite and errno == EDOM then ValueError is
995 raised
996 - if the result is finite and nonzero and errno == ERANGE then
997 OverflowError is raised
998
999 The last rule is used to catch overflow on platforms which follow
1000 C89 but for which HUGE_VAL is not an infinity.
1001
1002 For most two-argument functions (copysign, fmod, hypot, atan2)
1003 these rules are enough to ensure that Python's functions behave as
1004 specified in 'Annex F' of the C99 standard, with the 'invalid' and
1005 'divide-by-zero' floating-point exceptions mapping to Python's
1006 ValueError and the 'overflow' floating-point exception mapping to
1007 OverflowError.
1008 */
1009
1010 static PyObject *
math_1(PyObject * arg,double (* func)(double),int can_overflow)1011 math_1(PyObject *arg, double (*func) (double), int can_overflow)
1012 {
1013 return math_1_to_whatever(arg, func, PyFloat_FromDouble, can_overflow);
1014 }
1015
1016 static PyObject *
math_1_to_int(PyObject * arg,double (* func)(double),int can_overflow)1017 math_1_to_int(PyObject *arg, double (*func) (double), int can_overflow)
1018 {
1019 return math_1_to_whatever(arg, func, PyLong_FromDouble, can_overflow);
1020 }
1021
1022 static PyObject *
math_2(PyObject * const * args,Py_ssize_t nargs,double (* func)(double,double),const char * funcname)1023 math_2(PyObject *const *args, Py_ssize_t nargs,
1024 double (*func) (double, double), const char *funcname)
1025 {
1026 double x, y, r;
1027 if (!_PyArg_CheckPositional(funcname, nargs, 2, 2))
1028 return NULL;
1029 x = PyFloat_AsDouble(args[0]);
1030 if (x == -1.0 && PyErr_Occurred()) {
1031 return NULL;
1032 }
1033 y = PyFloat_AsDouble(args[1]);
1034 if (y == -1.0 && PyErr_Occurred()) {
1035 return NULL;
1036 }
1037 errno = 0;
1038 PyFPE_START_PROTECT("in math_2", return 0);
1039 r = (*func)(x, y);
1040 PyFPE_END_PROTECT(r);
1041 if (Py_IS_NAN(r)) {
1042 if (!Py_IS_NAN(x) && !Py_IS_NAN(y))
1043 errno = EDOM;
1044 else
1045 errno = 0;
1046 }
1047 else if (Py_IS_INFINITY(r)) {
1048 if (Py_IS_FINITE(x) && Py_IS_FINITE(y))
1049 errno = ERANGE;
1050 else
1051 errno = 0;
1052 }
1053 if (errno && is_error(r))
1054 return NULL;
1055 else
1056 return PyFloat_FromDouble(r);
1057 }
1058
1059 #define FUNC1(funcname, func, can_overflow, docstring) \
1060 static PyObject * math_##funcname(PyObject *self, PyObject *args) { \
1061 return math_1(args, func, can_overflow); \
1062 }\
1063 PyDoc_STRVAR(math_##funcname##_doc, docstring);
1064
1065 #define FUNC1A(funcname, func, docstring) \
1066 static PyObject * math_##funcname(PyObject *self, PyObject *args) { \
1067 return math_1a(args, func); \
1068 }\
1069 PyDoc_STRVAR(math_##funcname##_doc, docstring);
1070
1071 #define FUNC2(funcname, func, docstring) \
1072 static PyObject * math_##funcname(PyObject *self, PyObject *const *args, Py_ssize_t nargs) { \
1073 return math_2(args, nargs, func, #funcname); \
1074 }\
1075 PyDoc_STRVAR(math_##funcname##_doc, docstring);
1076
1077 FUNC1(acos, acos, 0,
1078 "acos($module, x, /)\n--\n\n"
1079 "Return the arc cosine (measured in radians) of x.")
1080 FUNC1(acosh, m_acosh, 0,
1081 "acosh($module, x, /)\n--\n\n"
1082 "Return the inverse hyperbolic cosine of x.")
1083 FUNC1(asin, asin, 0,
1084 "asin($module, x, /)\n--\n\n"
1085 "Return the arc sine (measured in radians) of x.")
1086 FUNC1(asinh, m_asinh, 0,
1087 "asinh($module, x, /)\n--\n\n"
1088 "Return the inverse hyperbolic sine of x.")
1089 FUNC1(atan, atan, 0,
1090 "atan($module, x, /)\n--\n\n"
1091 "Return the arc tangent (measured in radians) of x.")
1092 FUNC2(atan2, m_atan2,
1093 "atan2($module, y, x, /)\n--\n\n"
1094 "Return the arc tangent (measured in radians) of y/x.\n\n"
1095 "Unlike atan(y/x), the signs of both x and y are considered.")
1096 FUNC1(atanh, m_atanh, 0,
1097 "atanh($module, x, /)\n--\n\n"
1098 "Return the inverse hyperbolic tangent of x.")
1099
1100 /*[clinic input]
1101 math.ceil
1102
1103 x as number: object
1104 /
1105
1106 Return the ceiling of x as an Integral.
1107
1108 This is the smallest integer >= x.
1109 [clinic start generated code]*/
1110
1111 static PyObject *
math_ceil(PyObject * module,PyObject * number)1112 math_ceil(PyObject *module, PyObject *number)
1113 /*[clinic end generated code: output=6c3b8a78bc201c67 input=2725352806399cab]*/
1114 {
1115 _Py_IDENTIFIER(__ceil__);
1116 PyObject *method, *result;
1117
1118 method = _PyObject_LookupSpecial(number, &PyId___ceil__);
1119 if (method == NULL) {
1120 if (PyErr_Occurred())
1121 return NULL;
1122 return math_1_to_int(number, ceil, 0);
1123 }
1124 result = _PyObject_CallNoArg(method);
1125 Py_DECREF(method);
1126 return result;
1127 }
1128
1129 FUNC2(copysign, copysign,
1130 "copysign($module, x, y, /)\n--\n\n"
1131 "Return a float with the magnitude (absolute value) of x but the sign of y.\n\n"
1132 "On platforms that support signed zeros, copysign(1.0, -0.0)\n"
1133 "returns -1.0.\n")
1134 FUNC1(cos, cos, 0,
1135 "cos($module, x, /)\n--\n\n"
1136 "Return the cosine of x (measured in radians).")
1137 FUNC1(cosh, cosh, 1,
1138 "cosh($module, x, /)\n--\n\n"
1139 "Return the hyperbolic cosine of x.")
1140 FUNC1A(erf, m_erf,
1141 "erf($module, x, /)\n--\n\n"
1142 "Error function at x.")
1143 FUNC1A(erfc, m_erfc,
1144 "erfc($module, x, /)\n--\n\n"
1145 "Complementary error function at x.")
1146 FUNC1(exp, exp, 1,
1147 "exp($module, x, /)\n--\n\n"
1148 "Return e raised to the power of x.")
1149 FUNC1(expm1, m_expm1, 1,
1150 "expm1($module, x, /)\n--\n\n"
1151 "Return exp(x)-1.\n\n"
1152 "This function avoids the loss of precision involved in the direct "
1153 "evaluation of exp(x)-1 for small x.")
1154 FUNC1(fabs, fabs, 0,
1155 "fabs($module, x, /)\n--\n\n"
1156 "Return the absolute value of the float x.")
1157
1158 /*[clinic input]
1159 math.floor
1160
1161 x as number: object
1162 /
1163
1164 Return the floor of x as an Integral.
1165
1166 This is the largest integer <= x.
1167 [clinic start generated code]*/
1168
1169 static PyObject *
math_floor(PyObject * module,PyObject * number)1170 math_floor(PyObject *module, PyObject *number)
1171 /*[clinic end generated code: output=c6a65c4884884b8a input=63af6b5d7ebcc3d6]*/
1172 {
1173 _Py_IDENTIFIER(__floor__);
1174 PyObject *method, *result;
1175
1176 method = _PyObject_LookupSpecial(number, &PyId___floor__);
1177 if (method == NULL) {
1178 if (PyErr_Occurred())
1179 return NULL;
1180 return math_1_to_int(number, floor, 0);
1181 }
1182 result = _PyObject_CallNoArg(method);
1183 Py_DECREF(method);
1184 return result;
1185 }
1186
1187 FUNC1A(gamma, m_tgamma,
1188 "gamma($module, x, /)\n--\n\n"
1189 "Gamma function at x.")
1190 FUNC1A(lgamma, m_lgamma,
1191 "lgamma($module, x, /)\n--\n\n"
1192 "Natural logarithm of absolute value of Gamma function at x.")
1193 FUNC1(log1p, m_log1p, 0,
1194 "log1p($module, x, /)\n--\n\n"
1195 "Return the natural logarithm of 1+x (base e).\n\n"
1196 "The result is computed in a way which is accurate for x near zero.")
1197 FUNC2(remainder, m_remainder,
1198 "remainder($module, x, y, /)\n--\n\n"
1199 "Difference between x and the closest integer multiple of y.\n\n"
1200 "Return x - n*y where n*y is the closest integer multiple of y.\n"
1201 "In the case where x is exactly halfway between two multiples of\n"
1202 "y, the nearest even value of n is used. The result is always exact.")
1203 FUNC1(sin, sin, 0,
1204 "sin($module, x, /)\n--\n\n"
1205 "Return the sine of x (measured in radians).")
1206 FUNC1(sinh, sinh, 1,
1207 "sinh($module, x, /)\n--\n\n"
1208 "Return the hyperbolic sine of x.")
1209 FUNC1(sqrt, sqrt, 0,
1210 "sqrt($module, x, /)\n--\n\n"
1211 "Return the square root of x.")
1212 FUNC1(tan, tan, 0,
1213 "tan($module, x, /)\n--\n\n"
1214 "Return the tangent of x (measured in radians).")
1215 FUNC1(tanh, tanh, 0,
1216 "tanh($module, x, /)\n--\n\n"
1217 "Return the hyperbolic tangent of x.")
1218
1219 /* Precision summation function as msum() by Raymond Hettinger in
1220 <http://aspn.activestate.com/ASPN/Cookbook/Python/Recipe/393090>,
1221 enhanced with the exact partials sum and roundoff from Mark
1222 Dickinson's post at <http://bugs.python.org/file10357/msum4.py>.
1223 See those links for more details, proofs and other references.
1224
1225 Note 1: IEEE 754R floating point semantics are assumed,
1226 but the current implementation does not re-establish special
1227 value semantics across iterations (i.e. handling -Inf + Inf).
1228
1229 Note 2: No provision is made for intermediate overflow handling;
1230 therefore, sum([1e+308, 1e-308, 1e+308]) returns 1e+308 while
1231 sum([1e+308, 1e+308, 1e-308]) raises an OverflowError due to the
1232 overflow of the first partial sum.
1233
1234 Note 3: The intermediate values lo, yr, and hi are declared volatile so
1235 aggressive compilers won't algebraically reduce lo to always be exactly 0.0.
1236 Also, the volatile declaration forces the values to be stored in memory as
1237 regular doubles instead of extended long precision (80-bit) values. This
1238 prevents double rounding because any addition or subtraction of two doubles
1239 can be resolved exactly into double-sized hi and lo values. As long as the
1240 hi value gets forced into a double before yr and lo are computed, the extra
1241 bits in downstream extended precision operations (x87 for example) will be
1242 exactly zero and therefore can be losslessly stored back into a double,
1243 thereby preventing double rounding.
1244
1245 Note 4: A similar implementation is in Modules/cmathmodule.c.
1246 Be sure to update both when making changes.
1247
1248 Note 5: The signature of math.fsum() differs from builtins.sum()
1249 because the start argument doesn't make sense in the context of
1250 accurate summation. Since the partials table is collapsed before
1251 returning a result, sum(seq2, start=sum(seq1)) may not equal the
1252 accurate result returned by sum(itertools.chain(seq1, seq2)).
1253 */
1254
1255 #define NUM_PARTIALS 32 /* initial partials array size, on stack */
1256
1257 /* Extend the partials array p[] by doubling its size. */
1258 static int /* non-zero on error */
_fsum_realloc(double ** p_ptr,Py_ssize_t n,double * ps,Py_ssize_t * m_ptr)1259 _fsum_realloc(double **p_ptr, Py_ssize_t n,
1260 double *ps, Py_ssize_t *m_ptr)
1261 {
1262 void *v = NULL;
1263 Py_ssize_t m = *m_ptr;
1264
1265 m += m; /* double */
1266 if (n < m && (size_t)m < ((size_t)PY_SSIZE_T_MAX / sizeof(double))) {
1267 double *p = *p_ptr;
1268 if (p == ps) {
1269 v = PyMem_Malloc(sizeof(double) * m);
1270 if (v != NULL)
1271 memcpy(v, ps, sizeof(double) * n);
1272 }
1273 else
1274 v = PyMem_Realloc(p, sizeof(double) * m);
1275 }
1276 if (v == NULL) { /* size overflow or no memory */
1277 PyErr_SetString(PyExc_MemoryError, "math.fsum partials");
1278 return 1;
1279 }
1280 *p_ptr = (double*) v;
1281 *m_ptr = m;
1282 return 0;
1283 }
1284
1285 /* Full precision summation of a sequence of floats.
1286
1287 def msum(iterable):
1288 partials = [] # sorted, non-overlapping partial sums
1289 for x in iterable:
1290 i = 0
1291 for y in partials:
1292 if abs(x) < abs(y):
1293 x, y = y, x
1294 hi = x + y
1295 lo = y - (hi - x)
1296 if lo:
1297 partials[i] = lo
1298 i += 1
1299 x = hi
1300 partials[i:] = [x]
1301 return sum_exact(partials)
1302
1303 Rounded x+y stored in hi with the roundoff stored in lo. Together hi+lo
1304 are exactly equal to x+y. The inner loop applies hi/lo summation to each
1305 partial so that the list of partial sums remains exact.
1306
1307 Sum_exact() adds the partial sums exactly and correctly rounds the final
1308 result (using the round-half-to-even rule). The items in partials remain
1309 non-zero, non-special, non-overlapping and strictly increasing in
1310 magnitude, but possibly not all having the same sign.
1311
1312 Depends on IEEE 754 arithmetic guarantees and half-even rounding.
1313 */
1314
1315 /*[clinic input]
1316 math.fsum
1317
1318 seq: object
1319 /
1320
1321 Return an accurate floating point sum of values in the iterable seq.
1322
1323 Assumes IEEE-754 floating point arithmetic.
1324 [clinic start generated code]*/
1325
1326 static PyObject *
math_fsum(PyObject * module,PyObject * seq)1327 math_fsum(PyObject *module, PyObject *seq)
1328 /*[clinic end generated code: output=ba5c672b87fe34fc input=c51b7d8caf6f6e82]*/
1329 {
1330 PyObject *item, *iter, *sum = NULL;
1331 Py_ssize_t i, j, n = 0, m = NUM_PARTIALS;
1332 double x, y, t, ps[NUM_PARTIALS], *p = ps;
1333 double xsave, special_sum = 0.0, inf_sum = 0.0;
1334 volatile double hi, yr, lo;
1335
1336 iter = PyObject_GetIter(seq);
1337 if (iter == NULL)
1338 return NULL;
1339
1340 PyFPE_START_PROTECT("fsum", Py_DECREF(iter); return NULL)
1341
1342 for(;;) { /* for x in iterable */
1343 assert(0 <= n && n <= m);
1344 assert((m == NUM_PARTIALS && p == ps) ||
1345 (m > NUM_PARTIALS && p != NULL));
1346
1347 item = PyIter_Next(iter);
1348 if (item == NULL) {
1349 if (PyErr_Occurred())
1350 goto _fsum_error;
1351 break;
1352 }
1353 ASSIGN_DOUBLE(x, item, error_with_item);
1354 Py_DECREF(item);
1355
1356 xsave = x;
1357 for (i = j = 0; j < n; j++) { /* for y in partials */
1358 y = p[j];
1359 if (fabs(x) < fabs(y)) {
1360 t = x; x = y; y = t;
1361 }
1362 hi = x + y;
1363 yr = hi - x;
1364 lo = y - yr;
1365 if (lo != 0.0)
1366 p[i++] = lo;
1367 x = hi;
1368 }
1369
1370 n = i; /* ps[i:] = [x] */
1371 if (x != 0.0) {
1372 if (! Py_IS_FINITE(x)) {
1373 /* a nonfinite x could arise either as
1374 a result of intermediate overflow, or
1375 as a result of a nan or inf in the
1376 summands */
1377 if (Py_IS_FINITE(xsave)) {
1378 PyErr_SetString(PyExc_OverflowError,
1379 "intermediate overflow in fsum");
1380 goto _fsum_error;
1381 }
1382 if (Py_IS_INFINITY(xsave))
1383 inf_sum += xsave;
1384 special_sum += xsave;
1385 /* reset partials */
1386 n = 0;
1387 }
1388 else if (n >= m && _fsum_realloc(&p, n, ps, &m))
1389 goto _fsum_error;
1390 else
1391 p[n++] = x;
1392 }
1393 }
1394
1395 if (special_sum != 0.0) {
1396 if (Py_IS_NAN(inf_sum))
1397 PyErr_SetString(PyExc_ValueError,
1398 "-inf + inf in fsum");
1399 else
1400 sum = PyFloat_FromDouble(special_sum);
1401 goto _fsum_error;
1402 }
1403
1404 hi = 0.0;
1405 if (n > 0) {
1406 hi = p[--n];
1407 /* sum_exact(ps, hi) from the top, stop when the sum becomes
1408 inexact. */
1409 while (n > 0) {
1410 x = hi;
1411 y = p[--n];
1412 assert(fabs(y) < fabs(x));
1413 hi = x + y;
1414 yr = hi - x;
1415 lo = y - yr;
1416 if (lo != 0.0)
1417 break;
1418 }
1419 /* Make half-even rounding work across multiple partials.
1420 Needed so that sum([1e-16, 1, 1e16]) will round-up the last
1421 digit to two instead of down to zero (the 1e-16 makes the 1
1422 slightly closer to two). With a potential 1 ULP rounding
1423 error fixed-up, math.fsum() can guarantee commutativity. */
1424 if (n > 0 && ((lo < 0.0 && p[n-1] < 0.0) ||
1425 (lo > 0.0 && p[n-1] > 0.0))) {
1426 y = lo * 2.0;
1427 x = hi + y;
1428 yr = x - hi;
1429 if (y == yr)
1430 hi = x;
1431 }
1432 }
1433 sum = PyFloat_FromDouble(hi);
1434
1435 _fsum_error:
1436 PyFPE_END_PROTECT(hi)
1437 Py_DECREF(iter);
1438 if (p != ps)
1439 PyMem_Free(p);
1440 return sum;
1441
1442 error_with_item:
1443 Py_DECREF(item);
1444 goto _fsum_error;
1445 }
1446
1447 #undef NUM_PARTIALS
1448
1449
1450 /* Return the smallest integer k such that n < 2**k, or 0 if n == 0.
1451 * Equivalent to floor(lg(x))+1. Also equivalent to: bitwidth_of_type -
1452 * count_leading_zero_bits(x)
1453 */
1454
1455 /* XXX: This routine does more or less the same thing as
1456 * bits_in_digit() in Objects/longobject.c. Someday it would be nice to
1457 * consolidate them. On BSD, there's a library function called fls()
1458 * that we could use, and GCC provides __builtin_clz().
1459 */
1460
1461 static unsigned long
bit_length(unsigned long n)1462 bit_length(unsigned long n)
1463 {
1464 unsigned long len = 0;
1465 while (n != 0) {
1466 ++len;
1467 n >>= 1;
1468 }
1469 return len;
1470 }
1471
1472 static unsigned long
count_set_bits(unsigned long n)1473 count_set_bits(unsigned long n)
1474 {
1475 unsigned long count = 0;
1476 while (n != 0) {
1477 ++count;
1478 n &= n - 1; /* clear least significant bit */
1479 }
1480 return count;
1481 }
1482
1483 /* Integer square root
1484
1485 Given a nonnegative integer `n`, we want to compute the largest integer
1486 `a` for which `a * a <= n`, or equivalently the integer part of the exact
1487 square root of `n`.
1488
1489 We use an adaptive-precision pure-integer version of Newton's iteration. Given
1490 a positive integer `n`, the algorithm produces at each iteration an integer
1491 approximation `a` to the square root of `n >> s` for some even integer `s`,
1492 with `s` decreasing as the iterations progress. On the final iteration, `s` is
1493 zero and we have an approximation to the square root of `n` itself.
1494
1495 At every step, the approximation `a` is strictly within 1.0 of the true square
1496 root, so we have
1497
1498 (a - 1)**2 < (n >> s) < (a + 1)**2
1499
1500 After the final iteration, a check-and-correct step is needed to determine
1501 whether `a` or `a - 1` gives the desired integer square root of `n`.
1502
1503 The algorithm is remarkable in its simplicity. There's no need for a
1504 per-iteration check-and-correct step, and termination is straightforward: the
1505 number of iterations is known in advance (it's exactly `floor(log2(log2(n)))`
1506 for `n > 1`). The only tricky part of the correctness proof is in establishing
1507 that the bound `(a - 1)**2 < (n >> s) < (a + 1)**2` is maintained from one
1508 iteration to the next. A sketch of the proof of this is given below.
1509
1510 In addition to the proof sketch, a formal, computer-verified proof
1511 of correctness (using Lean) of an equivalent recursive algorithm can be found
1512 here:
1513
1514 https://github.com/mdickinson/snippets/blob/master/proofs/isqrt/src/isqrt.lean
1515
1516
1517 Here's Python code equivalent to the C implementation below:
1518
1519 def isqrt(n):
1520 """
1521 Return the integer part of the square root of the input.
1522 """
1523 n = operator.index(n)
1524
1525 if n < 0:
1526 raise ValueError("isqrt() argument must be nonnegative")
1527 if n == 0:
1528 return 0
1529
1530 c = (n.bit_length() - 1) // 2
1531 a = 1
1532 d = 0
1533 for s in reversed(range(c.bit_length())):
1534 # Loop invariant: (a-1)**2 < (n >> 2*(c - d)) < (a+1)**2
1535 e = d
1536 d = c >> s
1537 a = (a << d - e - 1) + (n >> 2*c - e - d + 1) // a
1538
1539 return a - (a*a > n)
1540
1541
1542 Sketch of proof of correctness
1543 ------------------------------
1544
1545 The delicate part of the correctness proof is showing that the loop invariant
1546 is preserved from one iteration to the next. That is, just before the line
1547
1548 a = (a << d - e - 1) + (n >> 2*c - e - d + 1) // a
1549
1550 is executed in the above code, we know that
1551
1552 (1) (a - 1)**2 < (n >> 2*(c - e)) < (a + 1)**2.
1553
1554 (since `e` is always the value of `d` from the previous iteration). We must
1555 prove that after that line is executed, we have
1556
1557 (a - 1)**2 < (n >> 2*(c - d)) < (a + 1)**2
1558
1559 To facilitate the proof, we make some changes of notation. Write `m` for
1560 `n >> 2*(c-d)`, and write `b` for the new value of `a`, so
1561
1562 b = (a << d - e - 1) + (n >> 2*c - e - d + 1) // a
1563
1564 or equivalently:
1565
1566 (2) b = (a << d - e - 1) + (m >> d - e + 1) // a
1567
1568 Then we can rewrite (1) as:
1569
1570 (3) (a - 1)**2 < (m >> 2*(d - e)) < (a + 1)**2
1571
1572 and we must show that (b - 1)**2 < m < (b + 1)**2.
1573
1574 From this point on, we switch to mathematical notation, so `/` means exact
1575 division rather than integer division and `^` is used for exponentiation. We
1576 use the `√` symbol for the exact square root. In (3), we can remove the
1577 implicit floor operation to give:
1578
1579 (4) (a - 1)^2 < m / 4^(d - e) < (a + 1)^2
1580
1581 Taking square roots throughout (4), scaling by `2^(d-e)`, and rearranging gives
1582
1583 (5) 0 <= | 2^(d-e)a - √m | < 2^(d-e)
1584
1585 Squaring and dividing through by `2^(d-e+1) a` gives
1586
1587 (6) 0 <= 2^(d-e-1) a + m / (2^(d-e+1) a) - √m < 2^(d-e-1) / a
1588
1589 We'll show below that `2^(d-e-1) <= a`. Given that, we can replace the
1590 right-hand side of (6) with `1`, and now replacing the central
1591 term `m / (2^(d-e+1) a)` with its floor in (6) gives
1592
1593 (7) -1 < 2^(d-e-1) a + m // 2^(d-e+1) a - √m < 1
1594
1595 Or equivalently, from (2):
1596
1597 (7) -1 < b - √m < 1
1598
1599 and rearranging gives that `(b-1)^2 < m < (b+1)^2`, which is what we needed
1600 to prove.
1601
1602 We're not quite done: we still have to prove the inequality `2^(d - e - 1) <=
1603 a` that was used to get line (7) above. From the definition of `c`, we have
1604 `4^c <= n`, which implies
1605
1606 (8) 4^d <= m
1607
1608 also, since `e == d >> 1`, `d` is at most `2e + 1`, from which it follows
1609 that `2d - 2e - 1 <= d` and hence that
1610
1611 (9) 4^(2d - 2e - 1) <= m
1612
1613 Dividing both sides by `4^(d - e)` gives
1614
1615 (10) 4^(d - e - 1) <= m / 4^(d - e)
1616
1617 But we know from (4) that `m / 4^(d-e) < (a + 1)^2`, hence
1618
1619 (11) 4^(d - e - 1) < (a + 1)^2
1620
1621 Now taking square roots of both sides and observing that both `2^(d-e-1)` and
1622 `a` are integers gives `2^(d - e - 1) <= a`, which is what we needed. This
1623 completes the proof sketch.
1624
1625 */
1626
1627
1628 /* Approximate square root of a large 64-bit integer.
1629
1630 Given `n` satisfying `2**62 <= n < 2**64`, return `a`
1631 satisfying `(a - 1)**2 < n < (a + 1)**2`. */
1632
1633 static uint64_t
_approximate_isqrt(uint64_t n)1634 _approximate_isqrt(uint64_t n)
1635 {
1636 uint32_t u = 1U + (n >> 62);
1637 u = (u << 1) + (n >> 59) / u;
1638 u = (u << 3) + (n >> 53) / u;
1639 u = (u << 7) + (n >> 41) / u;
1640 return (u << 15) + (n >> 17) / u;
1641 }
1642
1643 /*[clinic input]
1644 math.isqrt
1645
1646 n: object
1647 /
1648
1649 Return the integer part of the square root of the input.
1650 [clinic start generated code]*/
1651
1652 static PyObject *
math_isqrt(PyObject * module,PyObject * n)1653 math_isqrt(PyObject *module, PyObject *n)
1654 /*[clinic end generated code: output=35a6f7f980beab26 input=5b6e7ae4fa6c43d6]*/
1655 {
1656 int a_too_large, c_bit_length;
1657 size_t c, d;
1658 uint64_t m, u;
1659 PyObject *a = NULL, *b;
1660
1661 n = PyNumber_Index(n);
1662 if (n == NULL) {
1663 return NULL;
1664 }
1665
1666 if (_PyLong_Sign(n) < 0) {
1667 PyErr_SetString(
1668 PyExc_ValueError,
1669 "isqrt() argument must be nonnegative");
1670 goto error;
1671 }
1672 if (_PyLong_Sign(n) == 0) {
1673 Py_DECREF(n);
1674 return PyLong_FromLong(0);
1675 }
1676
1677 /* c = (n.bit_length() - 1) // 2 */
1678 c = _PyLong_NumBits(n);
1679 if (c == (size_t)(-1)) {
1680 goto error;
1681 }
1682 c = (c - 1U) / 2U;
1683
1684 /* Fast path: if c <= 31 then n < 2**64 and we can compute directly with a
1685 fast, almost branch-free algorithm. In the final correction, we use `u*u
1686 - 1 >= m` instead of the simpler `u*u > m` in order to get the correct
1687 result in the corner case where `u=2**32`. */
1688 if (c <= 31U) {
1689 m = (uint64_t)PyLong_AsUnsignedLongLong(n);
1690 Py_DECREF(n);
1691 if (m == (uint64_t)(-1) && PyErr_Occurred()) {
1692 return NULL;
1693 }
1694 u = _approximate_isqrt(m << (62U - 2U*c)) >> (31U - c);
1695 u -= u * u - 1U >= m;
1696 return PyLong_FromUnsignedLongLong((unsigned long long)u);
1697 }
1698
1699 /* Slow path: n >= 2**64. We perform the first five iterations in C integer
1700 arithmetic, then switch to using Python long integers. */
1701
1702 /* From n >= 2**64 it follows that c.bit_length() >= 6. */
1703 c_bit_length = 6;
1704 while ((c >> c_bit_length) > 0U) {
1705 ++c_bit_length;
1706 }
1707
1708 /* Initialise d and a. */
1709 d = c >> (c_bit_length - 5);
1710 b = _PyLong_Rshift(n, 2U*c - 62U);
1711 if (b == NULL) {
1712 goto error;
1713 }
1714 m = (uint64_t)PyLong_AsUnsignedLongLong(b);
1715 Py_DECREF(b);
1716 if (m == (uint64_t)(-1) && PyErr_Occurred()) {
1717 goto error;
1718 }
1719 u = _approximate_isqrt(m) >> (31U - d);
1720 a = PyLong_FromUnsignedLongLong((unsigned long long)u);
1721 if (a == NULL) {
1722 goto error;
1723 }
1724
1725 for (int s = c_bit_length - 6; s >= 0; --s) {
1726 PyObject *q;
1727 size_t e = d;
1728
1729 d = c >> s;
1730
1731 /* q = (n >> 2*c - e - d + 1) // a */
1732 q = _PyLong_Rshift(n, 2U*c - d - e + 1U);
1733 if (q == NULL) {
1734 goto error;
1735 }
1736 Py_SETREF(q, PyNumber_FloorDivide(q, a));
1737 if (q == NULL) {
1738 goto error;
1739 }
1740
1741 /* a = (a << d - 1 - e) + q */
1742 Py_SETREF(a, _PyLong_Lshift(a, d - 1U - e));
1743 if (a == NULL) {
1744 Py_DECREF(q);
1745 goto error;
1746 }
1747 Py_SETREF(a, PyNumber_Add(a, q));
1748 Py_DECREF(q);
1749 if (a == NULL) {
1750 goto error;
1751 }
1752 }
1753
1754 /* The correct result is either a or a - 1. Figure out which, and
1755 decrement a if necessary. */
1756
1757 /* a_too_large = n < a * a */
1758 b = PyNumber_Multiply(a, a);
1759 if (b == NULL) {
1760 goto error;
1761 }
1762 a_too_large = PyObject_RichCompareBool(n, b, Py_LT);
1763 Py_DECREF(b);
1764 if (a_too_large == -1) {
1765 goto error;
1766 }
1767
1768 if (a_too_large) {
1769 Py_SETREF(a, PyNumber_Subtract(a, _PyLong_One));
1770 }
1771 Py_DECREF(n);
1772 return a;
1773
1774 error:
1775 Py_XDECREF(a);
1776 Py_DECREF(n);
1777 return NULL;
1778 }
1779
1780 /* Divide-and-conquer factorial algorithm
1781 *
1782 * Based on the formula and pseudo-code provided at:
1783 * http://www.luschny.de/math/factorial/binarysplitfact.html
1784 *
1785 * Faster algorithms exist, but they're more complicated and depend on
1786 * a fast prime factorization algorithm.
1787 *
1788 * Notes on the algorithm
1789 * ----------------------
1790 *
1791 * factorial(n) is written in the form 2**k * m, with m odd. k and m are
1792 * computed separately, and then combined using a left shift.
1793 *
1794 * The function factorial_odd_part computes the odd part m (i.e., the greatest
1795 * odd divisor) of factorial(n), using the formula:
1796 *
1797 * factorial_odd_part(n) =
1798 *
1799 * product_{i >= 0} product_{0 < j <= n / 2**i, j odd} j
1800 *
1801 * Example: factorial_odd_part(20) =
1802 *
1803 * (1) *
1804 * (1) *
1805 * (1 * 3 * 5) *
1806 * (1 * 3 * 5 * 7 * 9)
1807 * (1 * 3 * 5 * 7 * 9 * 11 * 13 * 15 * 17 * 19)
1808 *
1809 * Here i goes from large to small: the first term corresponds to i=4 (any
1810 * larger i gives an empty product), and the last term corresponds to i=0.
1811 * Each term can be computed from the last by multiplying by the extra odd
1812 * numbers required: e.g., to get from the penultimate term to the last one,
1813 * we multiply by (11 * 13 * 15 * 17 * 19).
1814 *
1815 * To see a hint of why this formula works, here are the same numbers as above
1816 * but with the even parts (i.e., the appropriate powers of 2) included. For
1817 * each subterm in the product for i, we multiply that subterm by 2**i:
1818 *
1819 * factorial(20) =
1820 *
1821 * (16) *
1822 * (8) *
1823 * (4 * 12 * 20) *
1824 * (2 * 6 * 10 * 14 * 18) *
1825 * (1 * 3 * 5 * 7 * 9 * 11 * 13 * 15 * 17 * 19)
1826 *
1827 * The factorial_partial_product function computes the product of all odd j in
1828 * range(start, stop) for given start and stop. It's used to compute the
1829 * partial products like (11 * 13 * 15 * 17 * 19) in the example above. It
1830 * operates recursively, repeatedly splitting the range into two roughly equal
1831 * pieces until the subranges are small enough to be computed using only C
1832 * integer arithmetic.
1833 *
1834 * The two-valuation k (i.e., the exponent of the largest power of 2 dividing
1835 * the factorial) is computed independently in the main math_factorial
1836 * function. By standard results, its value is:
1837 *
1838 * two_valuation = n//2 + n//4 + n//8 + ....
1839 *
1840 * It can be shown (e.g., by complete induction on n) that two_valuation is
1841 * equal to n - count_set_bits(n), where count_set_bits(n) gives the number of
1842 * '1'-bits in the binary expansion of n.
1843 */
1844
1845 /* factorial_partial_product: Compute product(range(start, stop, 2)) using
1846 * divide and conquer. Assumes start and stop are odd and stop > start.
1847 * max_bits must be >= bit_length(stop - 2). */
1848
1849 static PyObject *
factorial_partial_product(unsigned long start,unsigned long stop,unsigned long max_bits)1850 factorial_partial_product(unsigned long start, unsigned long stop,
1851 unsigned long max_bits)
1852 {
1853 unsigned long midpoint, num_operands;
1854 PyObject *left = NULL, *right = NULL, *result = NULL;
1855
1856 /* If the return value will fit an unsigned long, then we can
1857 * multiply in a tight, fast loop where each multiply is O(1).
1858 * Compute an upper bound on the number of bits required to store
1859 * the answer.
1860 *
1861 * Storing some integer z requires floor(lg(z))+1 bits, which is
1862 * conveniently the value returned by bit_length(z). The
1863 * product x*y will require at most
1864 * bit_length(x) + bit_length(y) bits to store, based
1865 * on the idea that lg product = lg x + lg y.
1866 *
1867 * We know that stop - 2 is the largest number to be multiplied. From
1868 * there, we have: bit_length(answer) <= num_operands *
1869 * bit_length(stop - 2)
1870 */
1871
1872 num_operands = (stop - start) / 2;
1873 /* The "num_operands <= 8 * SIZEOF_LONG" check guards against the
1874 * unlikely case of an overflow in num_operands * max_bits. */
1875 if (num_operands <= 8 * SIZEOF_LONG &&
1876 num_operands * max_bits <= 8 * SIZEOF_LONG) {
1877 unsigned long j, total;
1878 for (total = start, j = start + 2; j < stop; j += 2)
1879 total *= j;
1880 return PyLong_FromUnsignedLong(total);
1881 }
1882
1883 /* find midpoint of range(start, stop), rounded up to next odd number. */
1884 midpoint = (start + num_operands) | 1;
1885 left = factorial_partial_product(start, midpoint,
1886 bit_length(midpoint - 2));
1887 if (left == NULL)
1888 goto error;
1889 right = factorial_partial_product(midpoint, stop, max_bits);
1890 if (right == NULL)
1891 goto error;
1892 result = PyNumber_Multiply(left, right);
1893
1894 error:
1895 Py_XDECREF(left);
1896 Py_XDECREF(right);
1897 return result;
1898 }
1899
1900 /* factorial_odd_part: compute the odd part of factorial(n). */
1901
1902 static PyObject *
factorial_odd_part(unsigned long n)1903 factorial_odd_part(unsigned long n)
1904 {
1905 long i;
1906 unsigned long v, lower, upper;
1907 PyObject *partial, *tmp, *inner, *outer;
1908
1909 inner = PyLong_FromLong(1);
1910 if (inner == NULL)
1911 return NULL;
1912 outer = inner;
1913 Py_INCREF(outer);
1914
1915 upper = 3;
1916 for (i = bit_length(n) - 2; i >= 0; i--) {
1917 v = n >> i;
1918 if (v <= 2)
1919 continue;
1920 lower = upper;
1921 /* (v + 1) | 1 = least odd integer strictly larger than n / 2**i */
1922 upper = (v + 1) | 1;
1923 /* Here inner is the product of all odd integers j in the range (0,
1924 n/2**(i+1)]. The factorial_partial_product call below gives the
1925 product of all odd integers j in the range (n/2**(i+1), n/2**i]. */
1926 partial = factorial_partial_product(lower, upper, bit_length(upper-2));
1927 /* inner *= partial */
1928 if (partial == NULL)
1929 goto error;
1930 tmp = PyNumber_Multiply(inner, partial);
1931 Py_DECREF(partial);
1932 if (tmp == NULL)
1933 goto error;
1934 Py_DECREF(inner);
1935 inner = tmp;
1936 /* Now inner is the product of all odd integers j in the range (0,
1937 n/2**i], giving the inner product in the formula above. */
1938
1939 /* outer *= inner; */
1940 tmp = PyNumber_Multiply(outer, inner);
1941 if (tmp == NULL)
1942 goto error;
1943 Py_DECREF(outer);
1944 outer = tmp;
1945 }
1946 Py_DECREF(inner);
1947 return outer;
1948
1949 error:
1950 Py_DECREF(outer);
1951 Py_DECREF(inner);
1952 return NULL;
1953 }
1954
1955
1956 /* Lookup table for small factorial values */
1957
1958 static const unsigned long SmallFactorials[] = {
1959 1, 1, 2, 6, 24, 120, 720, 5040, 40320,
1960 362880, 3628800, 39916800, 479001600,
1961 #if SIZEOF_LONG >= 8
1962 6227020800, 87178291200, 1307674368000,
1963 20922789888000, 355687428096000, 6402373705728000,
1964 121645100408832000, 2432902008176640000
1965 #endif
1966 };
1967
1968 /*[clinic input]
1969 math.factorial
1970
1971 x as arg: object
1972 /
1973
1974 Find x!.
1975
1976 Raise a ValueError if x is negative or non-integral.
1977 [clinic start generated code]*/
1978
1979 static PyObject *
math_factorial(PyObject * module,PyObject * arg)1980 math_factorial(PyObject *module, PyObject *arg)
1981 /*[clinic end generated code: output=6686f26fae00e9ca input=6d1c8105c0d91fb4]*/
1982 {
1983 long x, two_valuation;
1984 int overflow;
1985 PyObject *result, *odd_part, *pyint_form;
1986
1987 if (PyFloat_Check(arg)) {
1988 PyObject *lx;
1989 double dx = PyFloat_AS_DOUBLE((PyFloatObject *)arg);
1990 if (!(Py_IS_FINITE(dx) && dx == floor(dx))) {
1991 PyErr_SetString(PyExc_ValueError,
1992 "factorial() only accepts integral values");
1993 return NULL;
1994 }
1995 lx = PyLong_FromDouble(dx);
1996 if (lx == NULL)
1997 return NULL;
1998 x = PyLong_AsLongAndOverflow(lx, &overflow);
1999 Py_DECREF(lx);
2000 }
2001 else {
2002 pyint_form = PyNumber_Index(arg);
2003 if (pyint_form == NULL) {
2004 return NULL;
2005 }
2006 x = PyLong_AsLongAndOverflow(pyint_form, &overflow);
2007 Py_DECREF(pyint_form);
2008 }
2009
2010 if (x == -1 && PyErr_Occurred()) {
2011 return NULL;
2012 }
2013 else if (overflow == 1) {
2014 PyErr_Format(PyExc_OverflowError,
2015 "factorial() argument should not exceed %ld",
2016 LONG_MAX);
2017 return NULL;
2018 }
2019 else if (overflow == -1 || x < 0) {
2020 PyErr_SetString(PyExc_ValueError,
2021 "factorial() not defined for negative values");
2022 return NULL;
2023 }
2024
2025 /* use lookup table if x is small */
2026 if (x < (long)Py_ARRAY_LENGTH(SmallFactorials))
2027 return PyLong_FromUnsignedLong(SmallFactorials[x]);
2028
2029 /* else express in the form odd_part * 2**two_valuation, and compute as
2030 odd_part << two_valuation. */
2031 odd_part = factorial_odd_part(x);
2032 if (odd_part == NULL)
2033 return NULL;
2034 two_valuation = x - count_set_bits(x);
2035 result = _PyLong_Lshift(odd_part, two_valuation);
2036 Py_DECREF(odd_part);
2037 return result;
2038 }
2039
2040
2041 /*[clinic input]
2042 math.trunc
2043
2044 x: object
2045 /
2046
2047 Truncates the Real x to the nearest Integral toward 0.
2048
2049 Uses the __trunc__ magic method.
2050 [clinic start generated code]*/
2051
2052 static PyObject *
math_trunc(PyObject * module,PyObject * x)2053 math_trunc(PyObject *module, PyObject *x)
2054 /*[clinic end generated code: output=34b9697b707e1031 input=2168b34e0a09134d]*/
2055 {
2056 _Py_IDENTIFIER(__trunc__);
2057 PyObject *trunc, *result;
2058
2059 if (Py_TYPE(x)->tp_dict == NULL) {
2060 if (PyType_Ready(Py_TYPE(x)) < 0)
2061 return NULL;
2062 }
2063
2064 trunc = _PyObject_LookupSpecial(x, &PyId___trunc__);
2065 if (trunc == NULL) {
2066 if (!PyErr_Occurred())
2067 PyErr_Format(PyExc_TypeError,
2068 "type %.100s doesn't define __trunc__ method",
2069 Py_TYPE(x)->tp_name);
2070 return NULL;
2071 }
2072 result = _PyObject_CallNoArg(trunc);
2073 Py_DECREF(trunc);
2074 return result;
2075 }
2076
2077
2078 /*[clinic input]
2079 math.frexp
2080
2081 x: double
2082 /
2083
2084 Return the mantissa and exponent of x, as pair (m, e).
2085
2086 m is a float and e is an int, such that x = m * 2.**e.
2087 If x is 0, m and e are both 0. Else 0.5 <= abs(m) < 1.0.
2088 [clinic start generated code]*/
2089
2090 static PyObject *
math_frexp_impl(PyObject * module,double x)2091 math_frexp_impl(PyObject *module, double x)
2092 /*[clinic end generated code: output=03e30d252a15ad4a input=96251c9e208bc6e9]*/
2093 {
2094 int i;
2095 /* deal with special cases directly, to sidestep platform
2096 differences */
2097 if (Py_IS_NAN(x) || Py_IS_INFINITY(x) || !x) {
2098 i = 0;
2099 }
2100 else {
2101 PyFPE_START_PROTECT("in math_frexp", return 0);
2102 x = frexp(x, &i);
2103 PyFPE_END_PROTECT(x);
2104 }
2105 return Py_BuildValue("(di)", x, i);
2106 }
2107
2108
2109 /*[clinic input]
2110 math.ldexp
2111
2112 x: double
2113 i: object
2114 /
2115
2116 Return x * (2**i).
2117
2118 This is essentially the inverse of frexp().
2119 [clinic start generated code]*/
2120
2121 static PyObject *
math_ldexp_impl(PyObject * module,double x,PyObject * i)2122 math_ldexp_impl(PyObject *module, double x, PyObject *i)
2123 /*[clinic end generated code: output=b6892f3c2df9cc6a input=17d5970c1a40a8c1]*/
2124 {
2125 double r;
2126 long exp;
2127 int overflow;
2128
2129 if (PyLong_Check(i)) {
2130 /* on overflow, replace exponent with either LONG_MAX
2131 or LONG_MIN, depending on the sign. */
2132 exp = PyLong_AsLongAndOverflow(i, &overflow);
2133 if (exp == -1 && PyErr_Occurred())
2134 return NULL;
2135 if (overflow)
2136 exp = overflow < 0 ? LONG_MIN : LONG_MAX;
2137 }
2138 else {
2139 PyErr_SetString(PyExc_TypeError,
2140 "Expected an int as second argument to ldexp.");
2141 return NULL;
2142 }
2143
2144 if (x == 0. || !Py_IS_FINITE(x)) {
2145 /* NaNs, zeros and infinities are returned unchanged */
2146 r = x;
2147 errno = 0;
2148 } else if (exp > INT_MAX) {
2149 /* overflow */
2150 r = copysign(Py_HUGE_VAL, x);
2151 errno = ERANGE;
2152 } else if (exp < INT_MIN) {
2153 /* underflow to +-0 */
2154 r = copysign(0., x);
2155 errno = 0;
2156 } else {
2157 errno = 0;
2158 PyFPE_START_PROTECT("in math_ldexp", return 0);
2159 r = ldexp(x, (int)exp);
2160 PyFPE_END_PROTECT(r);
2161 if (Py_IS_INFINITY(r))
2162 errno = ERANGE;
2163 }
2164
2165 if (errno && is_error(r))
2166 return NULL;
2167 return PyFloat_FromDouble(r);
2168 }
2169
2170
2171 /*[clinic input]
2172 math.modf
2173
2174 x: double
2175 /
2176
2177 Return the fractional and integer parts of x.
2178
2179 Both results carry the sign of x and are floats.
2180 [clinic start generated code]*/
2181
2182 static PyObject *
math_modf_impl(PyObject * module,double x)2183 math_modf_impl(PyObject *module, double x)
2184 /*[clinic end generated code: output=90cee0260014c3c0 input=b4cfb6786afd9035]*/
2185 {
2186 double y;
2187 /* some platforms don't do the right thing for NaNs and
2188 infinities, so we take care of special cases directly. */
2189 if (!Py_IS_FINITE(x)) {
2190 if (Py_IS_INFINITY(x))
2191 return Py_BuildValue("(dd)", copysign(0., x), x);
2192 else if (Py_IS_NAN(x))
2193 return Py_BuildValue("(dd)", x, x);
2194 }
2195
2196 errno = 0;
2197 PyFPE_START_PROTECT("in math_modf", return 0);
2198 x = modf(x, &y);
2199 PyFPE_END_PROTECT(x);
2200 return Py_BuildValue("(dd)", x, y);
2201 }
2202
2203
2204 /* A decent logarithm is easy to compute even for huge ints, but libm can't
2205 do that by itself -- loghelper can. func is log or log10, and name is
2206 "log" or "log10". Note that overflow of the result isn't possible: an int
2207 can contain no more than INT_MAX * SHIFT bits, so has value certainly less
2208 than 2**(2**64 * 2**16) == 2**2**80, and log2 of that is 2**80, which is
2209 small enough to fit in an IEEE single. log and log10 are even smaller.
2210 However, intermediate overflow is possible for an int if the number of bits
2211 in that int is larger than PY_SSIZE_T_MAX. */
2212
2213 static PyObject*
loghelper(PyObject * arg,double (* func)(double),const char * funcname)2214 loghelper(PyObject* arg, double (*func)(double), const char *funcname)
2215 {
2216 /* If it is int, do it ourselves. */
2217 if (PyLong_Check(arg)) {
2218 double x, result;
2219 Py_ssize_t e;
2220
2221 /* Negative or zero inputs give a ValueError. */
2222 if (Py_SIZE(arg) <= 0) {
2223 PyErr_SetString(PyExc_ValueError,
2224 "math domain error");
2225 return NULL;
2226 }
2227
2228 x = PyLong_AsDouble(arg);
2229 if (x == -1.0 && PyErr_Occurred()) {
2230 if (!PyErr_ExceptionMatches(PyExc_OverflowError))
2231 return NULL;
2232 /* Here the conversion to double overflowed, but it's possible
2233 to compute the log anyway. Clear the exception and continue. */
2234 PyErr_Clear();
2235 x = _PyLong_Frexp((PyLongObject *)arg, &e);
2236 if (x == -1.0 && PyErr_Occurred())
2237 return NULL;
2238 /* Value is ~= x * 2**e, so the log ~= log(x) + log(2) * e. */
2239 result = func(x) + func(2.0) * e;
2240 }
2241 else
2242 /* Successfully converted x to a double. */
2243 result = func(x);
2244 return PyFloat_FromDouble(result);
2245 }
2246
2247 /* Else let libm handle it by itself. */
2248 return math_1(arg, func, 0);
2249 }
2250
2251
2252 /*[clinic input]
2253 math.log
2254
2255 x: object
2256 [
2257 base: object(c_default="NULL") = math.e
2258 ]
2259 /
2260
2261 Return the logarithm of x to the given base.
2262
2263 If the base not specified, returns the natural logarithm (base e) of x.
2264 [clinic start generated code]*/
2265
2266 static PyObject *
math_log_impl(PyObject * module,PyObject * x,int group_right_1,PyObject * base)2267 math_log_impl(PyObject *module, PyObject *x, int group_right_1,
2268 PyObject *base)
2269 /*[clinic end generated code: output=7b5a39e526b73fc9 input=0f62d5726cbfebbd]*/
2270 {
2271 PyObject *num, *den;
2272 PyObject *ans;
2273
2274 num = loghelper(x, m_log, "log");
2275 if (num == NULL || base == NULL)
2276 return num;
2277
2278 den = loghelper(base, m_log, "log");
2279 if (den == NULL) {
2280 Py_DECREF(num);
2281 return NULL;
2282 }
2283
2284 ans = PyNumber_TrueDivide(num, den);
2285 Py_DECREF(num);
2286 Py_DECREF(den);
2287 return ans;
2288 }
2289
2290
2291 /*[clinic input]
2292 math.log2
2293
2294 x: object
2295 /
2296
2297 Return the base 2 logarithm of x.
2298 [clinic start generated code]*/
2299
2300 static PyObject *
math_log2(PyObject * module,PyObject * x)2301 math_log2(PyObject *module, PyObject *x)
2302 /*[clinic end generated code: output=5425899a4d5d6acb input=08321262bae4f39b]*/
2303 {
2304 return loghelper(x, m_log2, "log2");
2305 }
2306
2307
2308 /*[clinic input]
2309 math.log10
2310
2311 x: object
2312 /
2313
2314 Return the base 10 logarithm of x.
2315 [clinic start generated code]*/
2316
2317 static PyObject *
math_log10(PyObject * module,PyObject * x)2318 math_log10(PyObject *module, PyObject *x)
2319 /*[clinic end generated code: output=be72a64617df9c6f input=b2469d02c6469e53]*/
2320 {
2321 return loghelper(x, m_log10, "log10");
2322 }
2323
2324
2325 /*[clinic input]
2326 math.fmod
2327
2328 x: double
2329 y: double
2330 /
2331
2332 Return fmod(x, y), according to platform C.
2333
2334 x % y may differ.
2335 [clinic start generated code]*/
2336
2337 static PyObject *
math_fmod_impl(PyObject * module,double x,double y)2338 math_fmod_impl(PyObject *module, double x, double y)
2339 /*[clinic end generated code: output=7559d794343a27b5 input=4f84caa8cfc26a03]*/
2340 {
2341 double r;
2342 /* fmod(x, +/-Inf) returns x for finite x. */
2343 if (Py_IS_INFINITY(y) && Py_IS_FINITE(x))
2344 return PyFloat_FromDouble(x);
2345 errno = 0;
2346 PyFPE_START_PROTECT("in math_fmod", return 0);
2347 r = fmod(x, y);
2348 PyFPE_END_PROTECT(r);
2349 if (Py_IS_NAN(r)) {
2350 if (!Py_IS_NAN(x) && !Py_IS_NAN(y))
2351 errno = EDOM;
2352 else
2353 errno = 0;
2354 }
2355 if (errno && is_error(r))
2356 return NULL;
2357 else
2358 return PyFloat_FromDouble(r);
2359 }
2360
2361 /*
2362 Given an *n* length *vec* of values and a value *max*, compute:
2363
2364 max * sqrt(sum((x / max) ** 2 for x in vec))
2365
2366 The value of the *max* variable must be non-negative and
2367 equal to the absolute value of the largest magnitude
2368 entry in the vector. If n==0, then *max* should be 0.0.
2369 If an infinity is present in the vec, *max* should be INF.
2370
2371 The *found_nan* variable indicates whether some member of
2372 the *vec* is a NaN.
2373
2374 To improve accuracy and to increase the number of cases where
2375 vector_norm() is commutative, we use a variant of Neumaier
2376 summation specialized to exploit that we always know that
2377 |csum| >= |x|.
2378
2379 The *csum* variable tracks the cumulative sum and *frac* tracks
2380 the cumulative fractional errors at each step. Since this
2381 variant assumes that |csum| >= |x| at each step, we establish
2382 the precondition by starting the accumulation from 1.0 which
2383 represents the largest possible value of (x/max)**2.
2384
2385 After the loop is finished, the initial 1.0 is subtracted out
2386 for a net zero effect on the final sum. Since *csum* will be
2387 greater than 1.0, the subtraction of 1.0 will not cause
2388 fractional digits to be dropped from *csum*.
2389
2390 */
2391
2392 static inline double
vector_norm(Py_ssize_t n,double * vec,double max,int found_nan)2393 vector_norm(Py_ssize_t n, double *vec, double max, int found_nan)
2394 {
2395 double x, csum = 1.0, oldcsum, frac = 0.0;
2396 Py_ssize_t i;
2397
2398 if (Py_IS_INFINITY(max)) {
2399 return max;
2400 }
2401 if (found_nan) {
2402 return Py_NAN;
2403 }
2404 if (max == 0.0 || n <= 1) {
2405 return max;
2406 }
2407 for (i=0 ; i < n ; i++) {
2408 x = vec[i];
2409 assert(Py_IS_FINITE(x) && fabs(x) <= max);
2410 x /= max;
2411 x = x*x;
2412 oldcsum = csum;
2413 csum += x;
2414 assert(csum >= x);
2415 frac += (oldcsum - csum) + x;
2416 }
2417 return max * sqrt(csum - 1.0 + frac);
2418 }
2419
2420 #define NUM_STACK_ELEMS 16
2421
2422 /*[clinic input]
2423 math.dist
2424
2425 p: object
2426 q: object
2427 /
2428
2429 Return the Euclidean distance between two points p and q.
2430
2431 The points should be specified as sequences (or iterables) of
2432 coordinates. Both inputs must have the same dimension.
2433
2434 Roughly equivalent to:
2435 sqrt(sum((px - qx) ** 2.0 for px, qx in zip(p, q)))
2436 [clinic start generated code]*/
2437
2438 static PyObject *
math_dist_impl(PyObject * module,PyObject * p,PyObject * q)2439 math_dist_impl(PyObject *module, PyObject *p, PyObject *q)
2440 /*[clinic end generated code: output=56bd9538d06bbcfe input=74e85e1b6092e68e]*/
2441 {
2442 PyObject *item;
2443 double max = 0.0;
2444 double x, px, qx, result;
2445 Py_ssize_t i, m, n;
2446 int found_nan = 0, p_allocated = 0, q_allocated = 0;
2447 double diffs_on_stack[NUM_STACK_ELEMS];
2448 double *diffs = diffs_on_stack;
2449
2450 if (!PyTuple_Check(p)) {
2451 p = PySequence_Tuple(p);
2452 if (p == NULL) {
2453 return NULL;
2454 }
2455 p_allocated = 1;
2456 }
2457 if (!PyTuple_Check(q)) {
2458 q = PySequence_Tuple(q);
2459 if (q == NULL) {
2460 if (p_allocated) {
2461 Py_DECREF(p);
2462 }
2463 return NULL;
2464 }
2465 q_allocated = 1;
2466 }
2467
2468 m = PyTuple_GET_SIZE(p);
2469 n = PyTuple_GET_SIZE(q);
2470 if (m != n) {
2471 PyErr_SetString(PyExc_ValueError,
2472 "both points must have the same number of dimensions");
2473 return NULL;
2474
2475 }
2476 if (n > NUM_STACK_ELEMS) {
2477 diffs = (double *) PyObject_Malloc(n * sizeof(double));
2478 if (diffs == NULL) {
2479 return PyErr_NoMemory();
2480 }
2481 }
2482 for (i=0 ; i<n ; i++) {
2483 item = PyTuple_GET_ITEM(p, i);
2484 ASSIGN_DOUBLE(px, item, error_exit);
2485 item = PyTuple_GET_ITEM(q, i);
2486 ASSIGN_DOUBLE(qx, item, error_exit);
2487 x = fabs(px - qx);
2488 diffs[i] = x;
2489 found_nan |= Py_IS_NAN(x);
2490 if (x > max) {
2491 max = x;
2492 }
2493 }
2494 result = vector_norm(n, diffs, max, found_nan);
2495 if (diffs != diffs_on_stack) {
2496 PyObject_Free(diffs);
2497 }
2498 if (p_allocated) {
2499 Py_DECREF(p);
2500 }
2501 if (q_allocated) {
2502 Py_DECREF(q);
2503 }
2504 return PyFloat_FromDouble(result);
2505
2506 error_exit:
2507 if (diffs != diffs_on_stack) {
2508 PyObject_Free(diffs);
2509 }
2510 if (p_allocated) {
2511 Py_DECREF(p);
2512 }
2513 if (q_allocated) {
2514 Py_DECREF(q);
2515 }
2516 return NULL;
2517 }
2518
2519 /* AC: cannot convert yet, waiting for *args support */
2520 static PyObject *
math_hypot(PyObject * self,PyObject * const * args,Py_ssize_t nargs)2521 math_hypot(PyObject *self, PyObject *const *args, Py_ssize_t nargs)
2522 {
2523 Py_ssize_t i;
2524 PyObject *item;
2525 double max = 0.0;
2526 double x, result;
2527 int found_nan = 0;
2528 double coord_on_stack[NUM_STACK_ELEMS];
2529 double *coordinates = coord_on_stack;
2530
2531 if (nargs > NUM_STACK_ELEMS) {
2532 coordinates = (double *) PyObject_Malloc(nargs * sizeof(double));
2533 if (coordinates == NULL) {
2534 return PyErr_NoMemory();
2535 }
2536 }
2537 for (i = 0; i < nargs; i++) {
2538 item = args[i];
2539 ASSIGN_DOUBLE(x, item, error_exit);
2540 x = fabs(x);
2541 coordinates[i] = x;
2542 found_nan |= Py_IS_NAN(x);
2543 if (x > max) {
2544 max = x;
2545 }
2546 }
2547 result = vector_norm(nargs, coordinates, max, found_nan);
2548 if (coordinates != coord_on_stack) {
2549 PyObject_Free(coordinates);
2550 }
2551 return PyFloat_FromDouble(result);
2552
2553 error_exit:
2554 if (coordinates != coord_on_stack) {
2555 PyObject_Free(coordinates);
2556 }
2557 return NULL;
2558 }
2559
2560 #undef NUM_STACK_ELEMS
2561
2562 PyDoc_STRVAR(math_hypot_doc,
2563 "hypot(*coordinates) -> value\n\n\
2564 Multidimensional Euclidean distance from the origin to a point.\n\
2565 \n\
2566 Roughly equivalent to:\n\
2567 sqrt(sum(x**2 for x in coordinates))\n\
2568 \n\
2569 For a two dimensional point (x, y), gives the hypotenuse\n\
2570 using the Pythagorean theorem: sqrt(x*x + y*y).\n\
2571 \n\
2572 For example, the hypotenuse of a 3/4/5 right triangle is:\n\
2573 \n\
2574 >>> hypot(3.0, 4.0)\n\
2575 5.0\n\
2576 ");
2577
2578 /* pow can't use math_2, but needs its own wrapper: the problem is
2579 that an infinite result can arise either as a result of overflow
2580 (in which case OverflowError should be raised) or as a result of
2581 e.g. 0.**-5. (for which ValueError needs to be raised.)
2582 */
2583
2584 /*[clinic input]
2585 math.pow
2586
2587 x: double
2588 y: double
2589 /
2590
2591 Return x**y (x to the power of y).
2592 [clinic start generated code]*/
2593
2594 static PyObject *
math_pow_impl(PyObject * module,double x,double y)2595 math_pow_impl(PyObject *module, double x, double y)
2596 /*[clinic end generated code: output=fff93e65abccd6b0 input=c26f1f6075088bfd]*/
2597 {
2598 double r;
2599 int odd_y;
2600
2601 /* deal directly with IEEE specials, to cope with problems on various
2602 platforms whose semantics don't exactly match C99 */
2603 r = 0.; /* silence compiler warning */
2604 if (!Py_IS_FINITE(x) || !Py_IS_FINITE(y)) {
2605 errno = 0;
2606 if (Py_IS_NAN(x))
2607 r = y == 0. ? 1. : x; /* NaN**0 = 1 */
2608 else if (Py_IS_NAN(y))
2609 r = x == 1. ? 1. : y; /* 1**NaN = 1 */
2610 else if (Py_IS_INFINITY(x)) {
2611 odd_y = Py_IS_FINITE(y) && fmod(fabs(y), 2.0) == 1.0;
2612 if (y > 0.)
2613 r = odd_y ? x : fabs(x);
2614 else if (y == 0.)
2615 r = 1.;
2616 else /* y < 0. */
2617 r = odd_y ? copysign(0., x) : 0.;
2618 }
2619 else if (Py_IS_INFINITY(y)) {
2620 if (fabs(x) == 1.0)
2621 r = 1.;
2622 else if (y > 0. && fabs(x) > 1.0)
2623 r = y;
2624 else if (y < 0. && fabs(x) < 1.0) {
2625 r = -y; /* result is +inf */
2626 if (x == 0.) /* 0**-inf: divide-by-zero */
2627 errno = EDOM;
2628 }
2629 else
2630 r = 0.;
2631 }
2632 }
2633 else {
2634 /* let libm handle finite**finite */
2635 errno = 0;
2636 PyFPE_START_PROTECT("in math_pow", return 0);
2637 r = pow(x, y);
2638 PyFPE_END_PROTECT(r);
2639 /* a NaN result should arise only from (-ve)**(finite
2640 non-integer); in this case we want to raise ValueError. */
2641 if (!Py_IS_FINITE(r)) {
2642 if (Py_IS_NAN(r)) {
2643 errno = EDOM;
2644 }
2645 /*
2646 an infinite result here arises either from:
2647 (A) (+/-0.)**negative (-> divide-by-zero)
2648 (B) overflow of x**y with x and y finite
2649 */
2650 else if (Py_IS_INFINITY(r)) {
2651 if (x == 0.)
2652 errno = EDOM;
2653 else
2654 errno = ERANGE;
2655 }
2656 }
2657 }
2658
2659 if (errno && is_error(r))
2660 return NULL;
2661 else
2662 return PyFloat_FromDouble(r);
2663 }
2664
2665
2666 static const double degToRad = Py_MATH_PI / 180.0;
2667 static const double radToDeg = 180.0 / Py_MATH_PI;
2668
2669 /*[clinic input]
2670 math.degrees
2671
2672 x: double
2673 /
2674
2675 Convert angle x from radians to degrees.
2676 [clinic start generated code]*/
2677
2678 static PyObject *
math_degrees_impl(PyObject * module,double x)2679 math_degrees_impl(PyObject *module, double x)
2680 /*[clinic end generated code: output=7fea78b294acd12f input=81e016555d6e3660]*/
2681 {
2682 return PyFloat_FromDouble(x * radToDeg);
2683 }
2684
2685
2686 /*[clinic input]
2687 math.radians
2688
2689 x: double
2690 /
2691
2692 Convert angle x from degrees to radians.
2693 [clinic start generated code]*/
2694
2695 static PyObject *
math_radians_impl(PyObject * module,double x)2696 math_radians_impl(PyObject *module, double x)
2697 /*[clinic end generated code: output=34daa47caf9b1590 input=91626fc489fe3d63]*/
2698 {
2699 return PyFloat_FromDouble(x * degToRad);
2700 }
2701
2702
2703 /*[clinic input]
2704 math.isfinite
2705
2706 x: double
2707 /
2708
2709 Return True if x is neither an infinity nor a NaN, and False otherwise.
2710 [clinic start generated code]*/
2711
2712 static PyObject *
math_isfinite_impl(PyObject * module,double x)2713 math_isfinite_impl(PyObject *module, double x)
2714 /*[clinic end generated code: output=8ba1f396440c9901 input=46967d254812e54a]*/
2715 {
2716 return PyBool_FromLong((long)Py_IS_FINITE(x));
2717 }
2718
2719
2720 /*[clinic input]
2721 math.isnan
2722
2723 x: double
2724 /
2725
2726 Return True if x is a NaN (not a number), and False otherwise.
2727 [clinic start generated code]*/
2728
2729 static PyObject *
math_isnan_impl(PyObject * module,double x)2730 math_isnan_impl(PyObject *module, double x)
2731 /*[clinic end generated code: output=f537b4d6df878c3e input=935891e66083f46a]*/
2732 {
2733 return PyBool_FromLong((long)Py_IS_NAN(x));
2734 }
2735
2736
2737 /*[clinic input]
2738 math.isinf
2739
2740 x: double
2741 /
2742
2743 Return True if x is a positive or negative infinity, and False otherwise.
2744 [clinic start generated code]*/
2745
2746 static PyObject *
math_isinf_impl(PyObject * module,double x)2747 math_isinf_impl(PyObject *module, double x)
2748 /*[clinic end generated code: output=9f00cbec4de7b06b input=32630e4212cf961f]*/
2749 {
2750 return PyBool_FromLong((long)Py_IS_INFINITY(x));
2751 }
2752
2753
2754 /*[clinic input]
2755 math.isclose -> bool
2756
2757 a: double
2758 b: double
2759 *
2760 rel_tol: double = 1e-09
2761 maximum difference for being considered "close", relative to the
2762 magnitude of the input values
2763 abs_tol: double = 0.0
2764 maximum difference for being considered "close", regardless of the
2765 magnitude of the input values
2766
2767 Determine whether two floating point numbers are close in value.
2768
2769 Return True if a is close in value to b, and False otherwise.
2770
2771 For the values to be considered close, the difference between them
2772 must be smaller than at least one of the tolerances.
2773
2774 -inf, inf and NaN behave similarly to the IEEE 754 Standard. That
2775 is, NaN is not close to anything, even itself. inf and -inf are
2776 only close to themselves.
2777 [clinic start generated code]*/
2778
2779 static int
math_isclose_impl(PyObject * module,double a,double b,double rel_tol,double abs_tol)2780 math_isclose_impl(PyObject *module, double a, double b, double rel_tol,
2781 double abs_tol)
2782 /*[clinic end generated code: output=b73070207511952d input=f28671871ea5bfba]*/
2783 {
2784 double diff = 0.0;
2785
2786 /* sanity check on the inputs */
2787 if (rel_tol < 0.0 || abs_tol < 0.0 ) {
2788 PyErr_SetString(PyExc_ValueError,
2789 "tolerances must be non-negative");
2790 return -1;
2791 }
2792
2793 if ( a == b ) {
2794 /* short circuit exact equality -- needed to catch two infinities of
2795 the same sign. And perhaps speeds things up a bit sometimes.
2796 */
2797 return 1;
2798 }
2799
2800 /* This catches the case of two infinities of opposite sign, or
2801 one infinity and one finite number. Two infinities of opposite
2802 sign would otherwise have an infinite relative tolerance.
2803 Two infinities of the same sign are caught by the equality check
2804 above.
2805 */
2806
2807 if (Py_IS_INFINITY(a) || Py_IS_INFINITY(b)) {
2808 return 0;
2809 }
2810
2811 /* now do the regular computation
2812 this is essentially the "weak" test from the Boost library
2813 */
2814
2815 diff = fabs(b - a);
2816
2817 return (((diff <= fabs(rel_tol * b)) ||
2818 (diff <= fabs(rel_tol * a))) ||
2819 (diff <= abs_tol));
2820 }
2821
2822 static inline int
_check_long_mult_overflow(long a,long b)2823 _check_long_mult_overflow(long a, long b) {
2824
2825 /* From Python2's int_mul code:
2826
2827 Integer overflow checking for * is painful: Python tried a couple ways, but
2828 they didn't work on all platforms, or failed in endcases (a product of
2829 -sys.maxint-1 has been a particular pain).
2830
2831 Here's another way:
2832
2833 The native long product x*y is either exactly right or *way* off, being
2834 just the last n bits of the true product, where n is the number of bits
2835 in a long (the delivered product is the true product plus i*2**n for
2836 some integer i).
2837
2838 The native double product (double)x * (double)y is subject to three
2839 rounding errors: on a sizeof(long)==8 box, each cast to double can lose
2840 info, and even on a sizeof(long)==4 box, the multiplication can lose info.
2841 But, unlike the native long product, it's not in *range* trouble: even
2842 if sizeof(long)==32 (256-bit longs), the product easily fits in the
2843 dynamic range of a double. So the leading 50 (or so) bits of the double
2844 product are correct.
2845
2846 We check these two ways against each other, and declare victory if they're
2847 approximately the same. Else, because the native long product is the only
2848 one that can lose catastrophic amounts of information, it's the native long
2849 product that must have overflowed.
2850
2851 */
2852
2853 long longprod = (long)((unsigned long)a * b);
2854 double doubleprod = (double)a * (double)b;
2855 double doubled_longprod = (double)longprod;
2856
2857 if (doubled_longprod == doubleprod) {
2858 return 0;
2859 }
2860
2861 const double diff = doubled_longprod - doubleprod;
2862 const double absdiff = diff >= 0.0 ? diff : -diff;
2863 const double absprod = doubleprod >= 0.0 ? doubleprod : -doubleprod;
2864
2865 if (32.0 * absdiff <= absprod) {
2866 return 0;
2867 }
2868
2869 return 1;
2870 }
2871
2872 /*[clinic input]
2873 math.prod
2874
2875 iterable: object
2876 /
2877 *
2878 start: object(c_default="NULL") = 1
2879
2880 Calculate the product of all the elements in the input iterable.
2881
2882 The default start value for the product is 1.
2883
2884 When the iterable is empty, return the start value. This function is
2885 intended specifically for use with numeric values and may reject
2886 non-numeric types.
2887 [clinic start generated code]*/
2888
2889 static PyObject *
math_prod_impl(PyObject * module,PyObject * iterable,PyObject * start)2890 math_prod_impl(PyObject *module, PyObject *iterable, PyObject *start)
2891 /*[clinic end generated code: output=36153bedac74a198 input=4c5ab0682782ed54]*/
2892 {
2893 PyObject *result = start;
2894 PyObject *temp, *item, *iter;
2895
2896 iter = PyObject_GetIter(iterable);
2897 if (iter == NULL) {
2898 return NULL;
2899 }
2900
2901 if (result == NULL) {
2902 result = PyLong_FromLong(1);
2903 if (result == NULL) {
2904 Py_DECREF(iter);
2905 return NULL;
2906 }
2907 } else {
2908 Py_INCREF(result);
2909 }
2910 #ifndef SLOW_PROD
2911 /* Fast paths for integers keeping temporary products in C.
2912 * Assumes all inputs are the same type.
2913 * If the assumption fails, default to use PyObjects instead.
2914 */
2915 if (PyLong_CheckExact(result)) {
2916 int overflow;
2917 long i_result = PyLong_AsLongAndOverflow(result, &overflow);
2918 /* If this already overflowed, don't even enter the loop. */
2919 if (overflow == 0) {
2920 Py_DECREF(result);
2921 result = NULL;
2922 }
2923 /* Loop over all the items in the iterable until we finish, we overflow
2924 * or we found a non integer element */
2925 while(result == NULL) {
2926 item = PyIter_Next(iter);
2927 if (item == NULL) {
2928 Py_DECREF(iter);
2929 if (PyErr_Occurred()) {
2930 return NULL;
2931 }
2932 return PyLong_FromLong(i_result);
2933 }
2934 if (PyLong_CheckExact(item)) {
2935 long b = PyLong_AsLongAndOverflow(item, &overflow);
2936 if (overflow == 0 && !_check_long_mult_overflow(i_result, b)) {
2937 long x = i_result * b;
2938 i_result = x;
2939 Py_DECREF(item);
2940 continue;
2941 }
2942 }
2943 /* Either overflowed or is not an int.
2944 * Restore real objects and process normally */
2945 result = PyLong_FromLong(i_result);
2946 if (result == NULL) {
2947 Py_DECREF(item);
2948 Py_DECREF(iter);
2949 return NULL;
2950 }
2951 temp = PyNumber_Multiply(result, item);
2952 Py_DECREF(result);
2953 Py_DECREF(item);
2954 result = temp;
2955 if (result == NULL) {
2956 Py_DECREF(iter);
2957 return NULL;
2958 }
2959 }
2960 }
2961
2962 /* Fast paths for floats keeping temporary products in C.
2963 * Assumes all inputs are the same type.
2964 * If the assumption fails, default to use PyObjects instead.
2965 */
2966 if (PyFloat_CheckExact(result)) {
2967 double f_result = PyFloat_AS_DOUBLE(result);
2968 Py_DECREF(result);
2969 result = NULL;
2970 while(result == NULL) {
2971 item = PyIter_Next(iter);
2972 if (item == NULL) {
2973 Py_DECREF(iter);
2974 if (PyErr_Occurred()) {
2975 return NULL;
2976 }
2977 return PyFloat_FromDouble(f_result);
2978 }
2979 if (PyFloat_CheckExact(item)) {
2980 f_result *= PyFloat_AS_DOUBLE(item);
2981 Py_DECREF(item);
2982 continue;
2983 }
2984 if (PyLong_CheckExact(item)) {
2985 long value;
2986 int overflow;
2987 value = PyLong_AsLongAndOverflow(item, &overflow);
2988 if (!overflow) {
2989 f_result *= (double)value;
2990 Py_DECREF(item);
2991 continue;
2992 }
2993 }
2994 result = PyFloat_FromDouble(f_result);
2995 if (result == NULL) {
2996 Py_DECREF(item);
2997 Py_DECREF(iter);
2998 return NULL;
2999 }
3000 temp = PyNumber_Multiply(result, item);
3001 Py_DECREF(result);
3002 Py_DECREF(item);
3003 result = temp;
3004 if (result == NULL) {
3005 Py_DECREF(iter);
3006 return NULL;
3007 }
3008 }
3009 }
3010 #endif
3011 /* Consume rest of the iterable (if any) that could not be handled
3012 * by specialized functions above.*/
3013 for(;;) {
3014 item = PyIter_Next(iter);
3015 if (item == NULL) {
3016 /* error, or end-of-sequence */
3017 if (PyErr_Occurred()) {
3018 Py_DECREF(result);
3019 result = NULL;
3020 }
3021 break;
3022 }
3023 temp = PyNumber_Multiply(result, item);
3024 Py_DECREF(result);
3025 Py_DECREF(item);
3026 result = temp;
3027 if (result == NULL)
3028 break;
3029 }
3030 Py_DECREF(iter);
3031 return result;
3032 }
3033
3034
3035 /*[clinic input]
3036 math.perm
3037
3038 n: object
3039 k: object = None
3040 /
3041
3042 Number of ways to choose k items from n items without repetition and with order.
3043
3044 Evaluates to n! / (n - k)! when k <= n and evaluates
3045 to zero when k > n.
3046
3047 If k is not specified or is None, then k defaults to n
3048 and the function returns n!.
3049
3050 Raises TypeError if either of the arguments are not integers.
3051 Raises ValueError if either of the arguments are negative.
3052 [clinic start generated code]*/
3053
3054 static PyObject *
math_perm_impl(PyObject * module,PyObject * n,PyObject * k)3055 math_perm_impl(PyObject *module, PyObject *n, PyObject *k)
3056 /*[clinic end generated code: output=e021a25469653e23 input=5311c5a00f359b53]*/
3057 {
3058 PyObject *result = NULL, *factor = NULL;
3059 int overflow, cmp;
3060 long long i, factors;
3061
3062 if (k == Py_None) {
3063 return math_factorial(module, n);
3064 }
3065 n = PyNumber_Index(n);
3066 if (n == NULL) {
3067 return NULL;
3068 }
3069 if (!PyLong_CheckExact(n)) {
3070 Py_SETREF(n, _PyLong_Copy((PyLongObject *)n));
3071 if (n == NULL) {
3072 return NULL;
3073 }
3074 }
3075 k = PyNumber_Index(k);
3076 if (k == NULL) {
3077 Py_DECREF(n);
3078 return NULL;
3079 }
3080 if (!PyLong_CheckExact(k)) {
3081 Py_SETREF(k, _PyLong_Copy((PyLongObject *)k));
3082 if (k == NULL) {
3083 Py_DECREF(n);
3084 return NULL;
3085 }
3086 }
3087
3088 if (Py_SIZE(n) < 0) {
3089 PyErr_SetString(PyExc_ValueError,
3090 "n must be a non-negative integer");
3091 goto error;
3092 }
3093 if (Py_SIZE(k) < 0) {
3094 PyErr_SetString(PyExc_ValueError,
3095 "k must be a non-negative integer");
3096 goto error;
3097 }
3098
3099 cmp = PyObject_RichCompareBool(n, k, Py_LT);
3100 if (cmp != 0) {
3101 if (cmp > 0) {
3102 result = PyLong_FromLong(0);
3103 goto done;
3104 }
3105 goto error;
3106 }
3107
3108 factors = PyLong_AsLongLongAndOverflow(k, &overflow);
3109 if (overflow > 0) {
3110 PyErr_Format(PyExc_OverflowError,
3111 "k must not exceed %lld",
3112 LLONG_MAX);
3113 goto error;
3114 }
3115 else if (factors == -1) {
3116 /* k is nonnegative, so a return value of -1 can only indicate error */
3117 goto error;
3118 }
3119
3120 if (factors == 0) {
3121 result = PyLong_FromLong(1);
3122 goto done;
3123 }
3124
3125 result = n;
3126 Py_INCREF(result);
3127 if (factors == 1) {
3128 goto done;
3129 }
3130
3131 factor = n;
3132 Py_INCREF(factor);
3133 for (i = 1; i < factors; ++i) {
3134 Py_SETREF(factor, PyNumber_Subtract(factor, _PyLong_One));
3135 if (factor == NULL) {
3136 goto error;
3137 }
3138 Py_SETREF(result, PyNumber_Multiply(result, factor));
3139 if (result == NULL) {
3140 goto error;
3141 }
3142 }
3143 Py_DECREF(factor);
3144
3145 done:
3146 Py_DECREF(n);
3147 Py_DECREF(k);
3148 return result;
3149
3150 error:
3151 Py_XDECREF(factor);
3152 Py_XDECREF(result);
3153 Py_DECREF(n);
3154 Py_DECREF(k);
3155 return NULL;
3156 }
3157
3158
3159 /*[clinic input]
3160 math.comb
3161
3162 n: object
3163 k: object
3164 /
3165
3166 Number of ways to choose k items from n items without repetition and without order.
3167
3168 Evaluates to n! / (k! * (n - k)!) when k <= n and evaluates
3169 to zero when k > n.
3170
3171 Also called the binomial coefficient because it is equivalent
3172 to the coefficient of k-th term in polynomial expansion of the
3173 expression (1 + x)**n.
3174
3175 Raises TypeError if either of the arguments are not integers.
3176 Raises ValueError if either of the arguments are negative.
3177
3178 [clinic start generated code]*/
3179
3180 static PyObject *
math_comb_impl(PyObject * module,PyObject * n,PyObject * k)3181 math_comb_impl(PyObject *module, PyObject *n, PyObject *k)
3182 /*[clinic end generated code: output=bd2cec8d854f3493 input=9a05315af2518709]*/
3183 {
3184 PyObject *result = NULL, *factor = NULL, *temp;
3185 int overflow, cmp;
3186 long long i, factors;
3187
3188 n = PyNumber_Index(n);
3189 if (n == NULL) {
3190 return NULL;
3191 }
3192 if (!PyLong_CheckExact(n)) {
3193 Py_SETREF(n, _PyLong_Copy((PyLongObject *)n));
3194 if (n == NULL) {
3195 return NULL;
3196 }
3197 }
3198 k = PyNumber_Index(k);
3199 if (k == NULL) {
3200 Py_DECREF(n);
3201 return NULL;
3202 }
3203 if (!PyLong_CheckExact(k)) {
3204 Py_SETREF(k, _PyLong_Copy((PyLongObject *)k));
3205 if (k == NULL) {
3206 Py_DECREF(n);
3207 return NULL;
3208 }
3209 }
3210
3211 if (Py_SIZE(n) < 0) {
3212 PyErr_SetString(PyExc_ValueError,
3213 "n must be a non-negative integer");
3214 goto error;
3215 }
3216 if (Py_SIZE(k) < 0) {
3217 PyErr_SetString(PyExc_ValueError,
3218 "k must be a non-negative integer");
3219 goto error;
3220 }
3221
3222 /* k = min(k, n - k) */
3223 temp = PyNumber_Subtract(n, k);
3224 if (temp == NULL) {
3225 goto error;
3226 }
3227 if (Py_SIZE(temp) < 0) {
3228 Py_DECREF(temp);
3229 result = PyLong_FromLong(0);
3230 goto done;
3231 }
3232 cmp = PyObject_RichCompareBool(temp, k, Py_LT);
3233 if (cmp > 0) {
3234 Py_SETREF(k, temp);
3235 }
3236 else {
3237 Py_DECREF(temp);
3238 if (cmp < 0) {
3239 goto error;
3240 }
3241 }
3242
3243 factors = PyLong_AsLongLongAndOverflow(k, &overflow);
3244 if (overflow > 0) {
3245 PyErr_Format(PyExc_OverflowError,
3246 "min(n - k, k) must not exceed %lld",
3247 LLONG_MAX);
3248 goto error;
3249 }
3250 if (factors == -1) {
3251 /* k is nonnegative, so a return value of -1 can only indicate error */
3252 goto error;
3253 }
3254
3255 if (factors == 0) {
3256 result = PyLong_FromLong(1);
3257 goto done;
3258 }
3259
3260 result = n;
3261 Py_INCREF(result);
3262 if (factors == 1) {
3263 goto done;
3264 }
3265
3266 factor = n;
3267 Py_INCREF(factor);
3268 for (i = 1; i < factors; ++i) {
3269 Py_SETREF(factor, PyNumber_Subtract(factor, _PyLong_One));
3270 if (factor == NULL) {
3271 goto error;
3272 }
3273 Py_SETREF(result, PyNumber_Multiply(result, factor));
3274 if (result == NULL) {
3275 goto error;
3276 }
3277
3278 temp = PyLong_FromUnsignedLongLong((unsigned long long)i + 1);
3279 if (temp == NULL) {
3280 goto error;
3281 }
3282 Py_SETREF(result, PyNumber_FloorDivide(result, temp));
3283 Py_DECREF(temp);
3284 if (result == NULL) {
3285 goto error;
3286 }
3287 }
3288 Py_DECREF(factor);
3289
3290 done:
3291 Py_DECREF(n);
3292 Py_DECREF(k);
3293 return result;
3294
3295 error:
3296 Py_XDECREF(factor);
3297 Py_XDECREF(result);
3298 Py_DECREF(n);
3299 Py_DECREF(k);
3300 return NULL;
3301 }
3302
3303
3304 static PyMethodDef math_methods[] = {
3305 {"acos", math_acos, METH_O, math_acos_doc},
3306 {"acosh", math_acosh, METH_O, math_acosh_doc},
3307 {"asin", math_asin, METH_O, math_asin_doc},
3308 {"asinh", math_asinh, METH_O, math_asinh_doc},
3309 {"atan", math_atan, METH_O, math_atan_doc},
3310 {"atan2", (PyCFunction)(void(*)(void))math_atan2, METH_FASTCALL, math_atan2_doc},
3311 {"atanh", math_atanh, METH_O, math_atanh_doc},
3312 MATH_CEIL_METHODDEF
3313 {"copysign", (PyCFunction)(void(*)(void))math_copysign, METH_FASTCALL, math_copysign_doc},
3314 {"cos", math_cos, METH_O, math_cos_doc},
3315 {"cosh", math_cosh, METH_O, math_cosh_doc},
3316 MATH_DEGREES_METHODDEF
3317 MATH_DIST_METHODDEF
3318 {"erf", math_erf, METH_O, math_erf_doc},
3319 {"erfc", math_erfc, METH_O, math_erfc_doc},
3320 {"exp", math_exp, METH_O, math_exp_doc},
3321 {"expm1", math_expm1, METH_O, math_expm1_doc},
3322 {"fabs", math_fabs, METH_O, math_fabs_doc},
3323 MATH_FACTORIAL_METHODDEF
3324 MATH_FLOOR_METHODDEF
3325 MATH_FMOD_METHODDEF
3326 MATH_FREXP_METHODDEF
3327 MATH_FSUM_METHODDEF
3328 {"gamma", math_gamma, METH_O, math_gamma_doc},
3329 MATH_GCD_METHODDEF
3330 {"hypot", (PyCFunction)(void(*)(void))math_hypot, METH_FASTCALL, math_hypot_doc},
3331 MATH_ISCLOSE_METHODDEF
3332 MATH_ISFINITE_METHODDEF
3333 MATH_ISINF_METHODDEF
3334 MATH_ISNAN_METHODDEF
3335 MATH_ISQRT_METHODDEF
3336 MATH_LDEXP_METHODDEF
3337 {"lgamma", math_lgamma, METH_O, math_lgamma_doc},
3338 MATH_LOG_METHODDEF
3339 {"log1p", math_log1p, METH_O, math_log1p_doc},
3340 MATH_LOG10_METHODDEF
3341 MATH_LOG2_METHODDEF
3342 MATH_MODF_METHODDEF
3343 MATH_POW_METHODDEF
3344 MATH_RADIANS_METHODDEF
3345 {"remainder", (PyCFunction)(void(*)(void))math_remainder, METH_FASTCALL, math_remainder_doc},
3346 {"sin", math_sin, METH_O, math_sin_doc},
3347 {"sinh", math_sinh, METH_O, math_sinh_doc},
3348 {"sqrt", math_sqrt, METH_O, math_sqrt_doc},
3349 {"tan", math_tan, METH_O, math_tan_doc},
3350 {"tanh", math_tanh, METH_O, math_tanh_doc},
3351 MATH_TRUNC_METHODDEF
3352 MATH_PROD_METHODDEF
3353 MATH_PERM_METHODDEF
3354 MATH_COMB_METHODDEF
3355 {NULL, NULL} /* sentinel */
3356 };
3357
3358
3359 PyDoc_STRVAR(module_doc,
3360 "This module provides access to the mathematical functions\n"
3361 "defined by the C standard.");
3362
3363
3364 static struct PyModuleDef mathmodule = {
3365 PyModuleDef_HEAD_INIT,
3366 "math",
3367 module_doc,
3368 -1,
3369 math_methods,
3370 NULL,
3371 NULL,
3372 NULL,
3373 NULL
3374 };
3375
3376 PyMODINIT_FUNC
PyInit_math(void)3377 PyInit_math(void)
3378 {
3379 PyObject *m;
3380
3381 m = PyModule_Create(&mathmodule);
3382 if (m == NULL)
3383 goto finally;
3384
3385 PyModule_AddObject(m, "pi", PyFloat_FromDouble(Py_MATH_PI));
3386 PyModule_AddObject(m, "e", PyFloat_FromDouble(Py_MATH_E));
3387 PyModule_AddObject(m, "tau", PyFloat_FromDouble(Py_MATH_TAU)); /* 2pi */
3388 PyModule_AddObject(m, "inf", PyFloat_FromDouble(m_inf()));
3389 #if !defined(PY_NO_SHORT_FLOAT_REPR) || defined(Py_NAN)
3390 PyModule_AddObject(m, "nan", PyFloat_FromDouble(m_nan()));
3391 #endif
3392
3393 finally:
3394 return m;
3395 }
3396