xref: /freebsd/lib/msun/src/k_tan.c (revision 81ad6265)
1 /* @(#)k_tan.c 1.5 04/04/22 SMI */
2 
3 /*
4  * ====================================================
5  * Copyright 2004 Sun Microsystems, Inc.  All Rights Reserved.
6  *
7  * Permission to use, copy, modify, and distribute this
8  * software is freely granted, provided that this notice
9  * is preserved.
10  * ====================================================
11  */
12 
13 /* INDENT OFF */
14 #include <sys/cdefs.h>
15 __FBSDID("$FreeBSD$");
16 
17 /* __kernel_tan( x, y, k )
18  * kernel tan function on ~[-pi/4, pi/4] (except on -0), pi/4 ~ 0.7854
19  * Input x is assumed to be bounded by ~pi/4 in magnitude.
20  * Input y is the tail of x.
21  * Input k indicates whether tan (if k = 1) or -1/tan (if k = -1) is returned.
22  *
23  * Algorithm
24  *	1. Since tan(-x) = -tan(x), we need only to consider positive x.
25  *	2. Callers must return tan(-0) = -0 without calling here since our
26  *	   odd polynomial is not evaluated in a way that preserves -0.
27  *	   Callers may do the optimization tan(x) ~ x for tiny x.
28  *	3. tan(x) is approximated by a odd polynomial of degree 27 on
29  *	   [0,0.67434]
30  *		  	         3             27
31  *	   	tan(x) ~ x + T1*x + ... + T13*x
32  *	   where
33  *
34  * 	        |tan(x)         2     4            26   |     -59.2
35  * 	        |----- - (1+T1*x +T2*x +.... +T13*x    )| <= 2
36  * 	        |  x 					|
37  *
38  *	   Note: tan(x+y) = tan(x) + tan'(x)*y
39  *		          ~ tan(x) + (1+x*x)*y
40  *	   Therefore, for better accuracy in computing tan(x+y), let
41  *		     3      2      2       2       2
42  *		r = x *(T2+x *(T3+x *(...+x *(T12+x *T13))))
43  *	   then
44  *		 		    3    2
45  *		tan(x+y) = x + (T1*x + (x *(r+y)+y))
46  *
47  *      4. For x in [0.67434,pi/4],  let y = pi/4 - x, then
48  *		tan(x) = tan(pi/4-y) = (1-tan(y))/(1+tan(y))
49  *		       = 1 - 2*(tan(y) - (tan(y)^2)/(1+tan(y)))
50  */
51 
52 #include "math.h"
53 #include "math_private.h"
54 static const double xxx[] = {
55 		 3.33333333333334091986e-01,	/* 3FD55555, 55555563 */
56 		 1.33333333333201242699e-01,	/* 3FC11111, 1110FE7A */
57 		 5.39682539762260521377e-02,	/* 3FABA1BA, 1BB341FE */
58 		 2.18694882948595424599e-02,	/* 3F9664F4, 8406D637 */
59 		 8.86323982359930005737e-03,	/* 3F8226E3, E96E8493 */
60 		 3.59207910759131235356e-03,	/* 3F6D6D22, C9560328 */
61 		 1.45620945432529025516e-03,	/* 3F57DBC8, FEE08315 */
62 		 5.88041240820264096874e-04,	/* 3F4344D8, F2F26501 */
63 		 2.46463134818469906812e-04,	/* 3F3026F7, 1A8D1068 */
64 		 7.81794442939557092300e-05,	/* 3F147E88, A03792A6 */
65 		 7.14072491382608190305e-05,	/* 3F12B80F, 32F0A7E9 */
66 		-1.85586374855275456654e-05,	/* BEF375CB, DB605373 */
67 		 2.59073051863633712884e-05,	/* 3EFB2A70, 74BF7AD4 */
68 /* one */	 1.00000000000000000000e+00,	/* 3FF00000, 00000000 */
69 /* pio4 */	 7.85398163397448278999e-01,	/* 3FE921FB, 54442D18 */
70 /* pio4lo */	 3.06161699786838301793e-17	/* 3C81A626, 33145C07 */
71 };
72 #define	one	xxx[13]
73 #define	pio4	xxx[14]
74 #define	pio4lo	xxx[15]
75 #define	T	xxx
76 /* INDENT ON */
77 
78 double
79 __kernel_tan(double x, double y, int iy) {
80 	double z, r, v, w, s;
81 	int32_t ix, hx;
82 
83 	GET_HIGH_WORD(hx,x);
84 	ix = hx & 0x7fffffff;			/* high word of |x| */
85 	if (ix >= 0x3FE59428) {	/* |x| >= 0.6744 */
86 		if (hx < 0) {
87 			x = -x;
88 			y = -y;
89 		}
90 		z = pio4 - x;
91 		w = pio4lo - y;
92 		x = z + w;
93 		y = 0.0;
94 	}
95 	z = x * x;
96 	w = z * z;
97 	/*
98 	 * Break x^5*(T[1]+x^2*T[2]+...) into
99 	 * x^5(T[1]+x^4*T[3]+...+x^20*T[11]) +
100 	 * x^5(x^2*(T[2]+x^4*T[4]+...+x^22*[T12]))
101 	 */
102 	r = T[1] + w * (T[3] + w * (T[5] + w * (T[7] + w * (T[9] +
103 		w * T[11]))));
104 	v = z * (T[2] + w * (T[4] + w * (T[6] + w * (T[8] + w * (T[10] +
105 		w * T[12])))));
106 	s = z * x;
107 	r = y + z * (s * (r + v) + y);
108 	r += T[0] * s;
109 	w = x + r;
110 	if (ix >= 0x3FE59428) {
111 		v = (double) iy;
112 		return (double) (1 - ((hx >> 30) & 2)) *
113 			(v - 2.0 * (x - (w * w / (w + v) - r)));
114 	}
115 	if (iy == 1)
116 		return w;
117 	else {
118 		/*
119 		 * if allow error up to 2 ulp, simply return
120 		 * -1.0 / (x+r) here
121 		 */
122 		/* compute -1.0 / (x+r) accurately */
123 		double a, t;
124 		z = w;
125 		SET_LOW_WORD(z,0);
126 		v = r - (z - x);	/* z+v = r+x */
127 		t = a = -1.0 / w;	/* a = -1.0/w */
128 		SET_LOW_WORD(t,0);
129 		s = 1.0 + t * z;
130 		return t + a * (s + t * v);
131 	}
132 }
133