xref: /freebsd/lib/msun/src/catrigf.c (revision 02e9120893770924227138ba49df1edb3896112a)
1 /*-
2  * SPDX-License-Identifier: BSD-2-Clause
3  *
4  * Copyright (c) 2012 Stephen Montgomery-Smith <stephen@FreeBSD.ORG>
5  * All rights reserved.
6  *
7  * Redistribution and use in source and binary forms, with or without
8  * modification, are permitted provided that the following conditions
9  * are met:
10  * 1. Redistributions of source code must retain the above copyright
11  *    notice, this list of conditions and the following disclaimer.
12  * 2. Redistributions in binary form must reproduce the above copyright
13  *    notice, this list of conditions and the following disclaimer in the
14  *    documentation and/or other materials provided with the distribution.
15  *
16  * THIS SOFTWARE IS PROVIDED BY THE AUTHOR AND CONTRIBUTORS ``AS IS'' AND
17  * ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
18  * IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE
19  * ARE DISCLAIMED.  IN NO EVENT SHALL THE AUTHOR OR CONTRIBUTORS BE LIABLE
20  * FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL
21  * DAMAGES (INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS
22  * OR SERVICES; LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION)
23  * HOWEVER CAUSED AND ON ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT
24  * LIABILITY, OR TORT (INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY
25  * OUT OF THE USE OF THIS SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF
26  * SUCH DAMAGE.
27  */
28 
29 /*
30  * The algorithm is very close to that in "Implementing the complex arcsine
31  * and arccosine functions using exception handling" by T. E. Hull, Thomas F.
32  * Fairgrieve, and Ping Tak Peter Tang, published in ACM Transactions on
33  * Mathematical Software, Volume 23 Issue 3, 1997, Pages 299-335,
34  * http://dl.acm.org/citation.cfm?id=275324.
35  *
36  * See catrig.c for complete comments.
37  *
38  * XXX comments were removed automatically, and even short ones on the right
39  * of statements were removed (all of them), contrary to normal style.  Only
40  * a few comments on the right of declarations remain.
41  */
42 
43 #include <sys/cdefs.h>
44 #include <complex.h>
45 #include <float.h>
46 
47 #include "math.h"
48 #include "math_private.h"
49 
50 #undef isinf
51 #define isinf(x)	(fabsf(x) == INFINITY)
52 #undef isnan
53 #define isnan(x)	((x) != (x))
54 #define	raise_inexact()	do { volatile float junk __unused = 1 + tiny; } while(0)
55 #undef signbit
56 #define signbit(x)	(__builtin_signbitf(x))
57 
58 static const float
59 A_crossover =		10,
60 B_crossover =		0.6417,
61 FOUR_SQRT_MIN =		0x1p-61,
62 QUARTER_SQRT_MAX =	0x1p61,
63 m_e =			2.7182818285e0,		/*  0xadf854.0p-22 */
64 m_ln2 =			6.9314718056e-1,	/*  0xb17218.0p-24 */
65 pio2_hi =		1.5707962513e0,		/*  0xc90fda.0p-23 */
66 RECIP_EPSILON =		1 / FLT_EPSILON,
67 SQRT_3_EPSILON =	5.9801995673e-4,	/*  0x9cc471.0p-34 */
68 SQRT_6_EPSILON =	8.4572793338e-4,	/*  0xddb3d7.0p-34 */
69 SQRT_MIN =		0x1p-63;
70 
71 static const volatile float
72 pio2_lo =		7.5497899549e-8,	/*  0xa22169.0p-47 */
73 tiny =			0x1p-100;
74 
75 static float complex clog_for_large_values(float complex z);
76 
77 static inline float
78 f(float a, float b, float hypot_a_b)
79 {
80 	if (b < 0)
81 		return ((hypot_a_b - b) / 2);
82 	if (b == 0)
83 		return (a / 2);
84 	return (a * a / (hypot_a_b + b) / 2);
85 }
86 
87 static inline void
88 do_hard_work(float x, float y, float *rx, int *B_is_usable, float *B,
89     float *sqrt_A2my2, float *new_y)
90 {
91 	float R, S, A;
92 	float Am1, Amy;
93 
94 	R = hypotf(x, y + 1);
95 	S = hypotf(x, y - 1);
96 
97 	A = (R + S) / 2;
98 	if (A < 1)
99 		A = 1;
100 
101 	if (A < A_crossover) {
102 		if (y == 1 && x < FLT_EPSILON * FLT_EPSILON / 128) {
103 			*rx = sqrtf(x);
104 		} else if (x >= FLT_EPSILON * fabsf(y - 1)) {
105 			Am1 = f(x, 1 + y, R) + f(x, 1 - y, S);
106 			*rx = log1pf(Am1 + sqrtf(Am1 * (A + 1)));
107 		} else if (y < 1) {
108 			*rx = x / sqrtf((1 - y) * (1 + y));
109 		} else {
110 			*rx = log1pf((y - 1) + sqrtf((y - 1) * (y + 1)));
111 		}
112 	} else {
113 		*rx = logf(A + sqrtf(A * A - 1));
114 	}
115 
116 	*new_y = y;
117 
118 	if (y < FOUR_SQRT_MIN) {
119 		*B_is_usable = 0;
120 		*sqrt_A2my2 = A * (2 / FLT_EPSILON);
121 		*new_y = y * (2 / FLT_EPSILON);
122 		return;
123 	}
124 
125 	*B = y / A;
126 	*B_is_usable = 1;
127 
128 	if (*B > B_crossover) {
129 		*B_is_usable = 0;
130 		if (y == 1 && x < FLT_EPSILON / 128) {
131 			*sqrt_A2my2 = sqrtf(x) * sqrtf((A + y) / 2);
132 		} else if (x >= FLT_EPSILON * fabsf(y - 1)) {
133 			Amy = f(x, y + 1, R) + f(x, y - 1, S);
134 			*sqrt_A2my2 = sqrtf(Amy * (A + y));
135 		} else if (y > 1) {
136 			*sqrt_A2my2 = x * (4 / FLT_EPSILON / FLT_EPSILON) * y /
137 			    sqrtf((y + 1) * (y - 1));
138 			*new_y = y * (4 / FLT_EPSILON / FLT_EPSILON);
139 		} else {
140 			*sqrt_A2my2 = sqrtf((1 - y) * (1 + y));
141 		}
142 	}
143 }
144 
145 float complex
146 casinhf(float complex z)
147 {
148 	float x, y, ax, ay, rx, ry, B, sqrt_A2my2, new_y;
149 	int B_is_usable;
150 	float complex w;
151 
152 	x = crealf(z);
153 	y = cimagf(z);
154 	ax = fabsf(x);
155 	ay = fabsf(y);
156 
157 	if (isnan(x) || isnan(y)) {
158 		if (isinf(x))
159 			return (CMPLXF(x, y + y));
160 		if (isinf(y))
161 			return (CMPLXF(y, x + x));
162 		if (y == 0)
163 			return (CMPLXF(x + x, y));
164 		return (CMPLXF(nan_mix(x, y), nan_mix(x, y)));
165 	}
166 
167 	if (ax > RECIP_EPSILON || ay > RECIP_EPSILON) {
168 		if (signbit(x) == 0)
169 			w = clog_for_large_values(z) + m_ln2;
170 		else
171 			w = clog_for_large_values(-z) + m_ln2;
172 		return (CMPLXF(copysignf(crealf(w), x),
173 		    copysignf(cimagf(w), y)));
174 	}
175 
176 	if (x == 0 && y == 0)
177 		return (z);
178 
179 	raise_inexact();
180 
181 	if (ax < SQRT_6_EPSILON / 4 && ay < SQRT_6_EPSILON / 4)
182 		return (z);
183 
184 	do_hard_work(ax, ay, &rx, &B_is_usable, &B, &sqrt_A2my2, &new_y);
185 	if (B_is_usable)
186 		ry = asinf(B);
187 	else
188 		ry = atan2f(new_y, sqrt_A2my2);
189 	return (CMPLXF(copysignf(rx, x), copysignf(ry, y)));
190 }
191 
192 float complex
193 casinf(float complex z)
194 {
195 	float complex w = casinhf(CMPLXF(cimagf(z), crealf(z)));
196 
197 	return (CMPLXF(cimagf(w), crealf(w)));
198 }
199 
200 float complex
201 cacosf(float complex z)
202 {
203 	float x, y, ax, ay, rx, ry, B, sqrt_A2mx2, new_x;
204 	int sx, sy;
205 	int B_is_usable;
206 	float complex w;
207 
208 	x = crealf(z);
209 	y = cimagf(z);
210 	sx = signbit(x);
211 	sy = signbit(y);
212 	ax = fabsf(x);
213 	ay = fabsf(y);
214 
215 	if (isnan(x) || isnan(y)) {
216 		if (isinf(x))
217 			return (CMPLXF(y + y, -INFINITY));
218 		if (isinf(y))
219 			return (CMPLXF(x + x, -y));
220 		if (x == 0)
221 			return (CMPLXF(pio2_hi + pio2_lo, y + y));
222 		return (CMPLXF(nan_mix(x, y), nan_mix(x, y)));
223 	}
224 
225 	if (ax > RECIP_EPSILON || ay > RECIP_EPSILON) {
226 		w = clog_for_large_values(z);
227 		rx = fabsf(cimagf(w));
228 		ry = crealf(w) + m_ln2;
229 		if (sy == 0)
230 			ry = -ry;
231 		return (CMPLXF(rx, ry));
232 	}
233 
234 	if (x == 1 && y == 0)
235 		return (CMPLXF(0, -y));
236 
237 	raise_inexact();
238 
239 	if (ax < SQRT_6_EPSILON / 4 && ay < SQRT_6_EPSILON / 4)
240 		return (CMPLXF(pio2_hi - (x - pio2_lo), -y));
241 
242 	do_hard_work(ay, ax, &ry, &B_is_usable, &B, &sqrt_A2mx2, &new_x);
243 	if (B_is_usable) {
244 		if (sx == 0)
245 			rx = acosf(B);
246 		else
247 			rx = acosf(-B);
248 	} else {
249 		if (sx == 0)
250 			rx = atan2f(sqrt_A2mx2, new_x);
251 		else
252 			rx = atan2f(sqrt_A2mx2, -new_x);
253 	}
254 	if (sy == 0)
255 		ry = -ry;
256 	return (CMPLXF(rx, ry));
257 }
258 
259 float complex
260 cacoshf(float complex z)
261 {
262 	float complex w;
263 	float rx, ry;
264 
265 	w = cacosf(z);
266 	rx = crealf(w);
267 	ry = cimagf(w);
268 	if (isnan(rx) && isnan(ry))
269 		return (CMPLXF(ry, rx));
270 	if (isnan(rx))
271 		return (CMPLXF(fabsf(ry), rx));
272 	if (isnan(ry))
273 		return (CMPLXF(ry, ry));
274 	return (CMPLXF(fabsf(ry), copysignf(rx, cimagf(z))));
275 }
276 
277 static float complex
278 clog_for_large_values(float complex z)
279 {
280 	float x, y;
281 	float ax, ay, t;
282 
283 	x = crealf(z);
284 	y = cimagf(z);
285 	ax = fabsf(x);
286 	ay = fabsf(y);
287 	if (ax < ay) {
288 		t = ax;
289 		ax = ay;
290 		ay = t;
291 	}
292 
293 	if (ax > FLT_MAX / 2)
294 		return (CMPLXF(logf(hypotf(x / m_e, y / m_e)) + 1,
295 		    atan2f(y, x)));
296 
297 	if (ax > QUARTER_SQRT_MAX || ay < SQRT_MIN)
298 		return (CMPLXF(logf(hypotf(x, y)), atan2f(y, x)));
299 
300 	return (CMPLXF(logf(ax * ax + ay * ay) / 2, atan2f(y, x)));
301 }
302 
303 static inline float
304 sum_squares(float x, float y)
305 {
306 
307 	if (y < SQRT_MIN)
308 		return (x * x);
309 
310 	return (x * x + y * y);
311 }
312 
313 static inline float
314 real_part_reciprocal(float x, float y)
315 {
316 	float scale;
317 	uint32_t hx, hy;
318 	int32_t ix, iy;
319 
320 	GET_FLOAT_WORD(hx, x);
321 	ix = hx & 0x7f800000;
322 	GET_FLOAT_WORD(hy, y);
323 	iy = hy & 0x7f800000;
324 #define	BIAS	(FLT_MAX_EXP - 1)
325 #define	CUTOFF	(FLT_MANT_DIG / 2 + 1)
326 	if (ix - iy >= CUTOFF << 23 || isinf(x))
327 		return (1 / x);
328 	if (iy - ix >= CUTOFF << 23)
329 		return (x / y / y);
330 	if (ix <= (BIAS + FLT_MAX_EXP / 2 - CUTOFF) << 23)
331 		return (x / (x * x + y * y));
332 	SET_FLOAT_WORD(scale, 0x7f800000 - ix);
333 	x *= scale;
334 	y *= scale;
335 	return (x / (x * x + y * y) * scale);
336 }
337 
338 float complex
339 catanhf(float complex z)
340 {
341 	float x, y, ax, ay, rx, ry;
342 
343 	x = crealf(z);
344 	y = cimagf(z);
345 	ax = fabsf(x);
346 	ay = fabsf(y);
347 
348 	if (y == 0 && ax <= 1)
349 		return (CMPLXF(atanhf(x), y));
350 
351 	if (x == 0)
352 		return (CMPLXF(x, atanf(y)));
353 
354 	if (isnan(x) || isnan(y)) {
355 		if (isinf(x))
356 			return (CMPLXF(copysignf(0, x), y + y));
357 		if (isinf(y))
358 			return (CMPLXF(copysignf(0, x),
359 			    copysignf(pio2_hi + pio2_lo, y)));
360 		return (CMPLXF(nan_mix(x, y), nan_mix(x, y)));
361 	}
362 
363 	if (ax > RECIP_EPSILON || ay > RECIP_EPSILON)
364 		return (CMPLXF(real_part_reciprocal(x, y),
365 		    copysignf(pio2_hi + pio2_lo, y)));
366 
367 	if (ax < SQRT_3_EPSILON / 2 && ay < SQRT_3_EPSILON / 2) {
368 		raise_inexact();
369 		return (z);
370 	}
371 
372 	if (ax == 1 && ay < FLT_EPSILON)
373 		rx = (m_ln2 - logf(ay)) / 2;
374 	else
375 		rx = log1pf(4 * ax / sum_squares(ax - 1, ay)) / 4;
376 
377 	if (ax == 1)
378 		ry = atan2f(2, -ay) / 2;
379 	else if (ay < FLT_EPSILON)
380 		ry = atan2f(2 * ay, (1 - ax) * (1 + ax)) / 2;
381 	else
382 		ry = atan2f(2 * ay, (1 - ax) * (1 + ax) - ay * ay) / 2;
383 
384 	return (CMPLXF(copysignf(rx, x), copysignf(ry, y)));
385 }
386 
387 float complex
388 catanf(float complex z)
389 {
390 	float complex w = catanhf(CMPLXF(cimagf(z), crealf(z)));
391 
392 	return (CMPLXF(cimagf(w), crealf(w)));
393 }
394