Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
GenericPacketMathFunctions.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2007 Julien Pommier
5// Copyright (C) 2014 Pedro Gonnet ([email protected])
6// Copyright (C) 2009-2019 Gael Guennebaud <[email protected]>
7//
8// This Source Code Form is subject to the terms of the Mozilla
9// Public License v. 2.0. If a copy of the MPL was not distributed
10// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
11
12/* The exp and log functions of this file initially come from
13 * Julien Pommier's sse math library: http://gruntthepeon.free.fr/ssemath/
14 */
15
16#ifndef EIGEN_ARCH_GENERIC_PACKET_MATH_FUNCTIONS_H
17#define EIGEN_ARCH_GENERIC_PACKET_MATH_FUNCTIONS_H
18
19// IWYU pragma: private
20#include "../../InternalHeaderCheck.h"
21
22namespace Eigen {
23namespace internal {
24
25// Creates a Scalar integer type with same bit-width.
26template <typename T>
27struct make_integer;
28template <>
29struct make_integer<float> {
30 typedef numext::int32_t type;
31};
32template <>
33struct make_integer<double> {
34 typedef numext::int64_t type;
35};
36template <>
37struct make_integer<half> {
38 typedef numext::int16_t type;
39};
40template <>
41struct make_integer<bfloat16> {
42 typedef numext::int16_t type;
43};
44
45template <typename Packet>
46EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC Packet pfrexp_generic_get_biased_exponent(const Packet& a) {
47 typedef typename unpacket_traits<Packet>::type Scalar;
48 typedef typename unpacket_traits<Packet>::integer_packet PacketI;
49 static constexpr int mantissa_bits = numext::numeric_limits<Scalar>::digits - 1;
50 return pcast<PacketI, Packet>(plogical_shift_right<mantissa_bits>(preinterpret<PacketI>(pabs(a))));
51}
52
53// Safely applies frexp, correctly handles denormals.
54// Assumes IEEE floating point format.
55template <typename Packet>
56EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC Packet pfrexp_generic(const Packet& a, Packet& exponent) {
57 typedef typename unpacket_traits<Packet>::type Scalar;
58 typedef typename make_unsigned<typename make_integer<Scalar>::type>::type ScalarUI;
59 static constexpr int TotalBits = sizeof(Scalar) * CHAR_BIT, MantissaBits = numext::numeric_limits<Scalar>::digits - 1,
60 ExponentBits = TotalBits - MantissaBits - 1;
61
62 EIGEN_CONSTEXPR ScalarUI scalar_sign_mantissa_mask =
63 ~(((ScalarUI(1) << ExponentBits) - ScalarUI(1)) << MantissaBits); // ~0x7f800000
64 const Packet sign_mantissa_mask = pset1frombits<Packet>(static_cast<ScalarUI>(scalar_sign_mantissa_mask));
65 const Packet half = pset1<Packet>(Scalar(0.5));
66 const Packet zero = pzero(a);
67 const Packet normal_min = pset1<Packet>((numext::numeric_limits<Scalar>::min)()); // Minimum normal value, 2^-126
68
69 // To handle denormals, normalize by multiplying by 2^(int(MantissaBits)+1).
70 const Packet is_denormal = pcmp_lt(pabs(a), normal_min);
71 EIGEN_CONSTEXPR ScalarUI scalar_normalization_offset = ScalarUI(MantissaBits + 1); // 24
72 // The following cannot be constexpr because bfloat16(uint16_t) is not constexpr.
73 const Scalar scalar_normalization_factor = Scalar(ScalarUI(1) << int(scalar_normalization_offset)); // 2^24
74 const Packet normalization_factor = pset1<Packet>(scalar_normalization_factor);
75 const Packet normalized_a = pselect(is_denormal, pmul(a, normalization_factor), a);
76
77 // Determine exponent offset: -126 if normal, -126-24 if denormal
78 const Scalar scalar_exponent_offset = -Scalar((ScalarUI(1) << (ExponentBits - 1)) - ScalarUI(2)); // -126
79 Packet exponent_offset = pset1<Packet>(scalar_exponent_offset);
80 const Packet normalization_offset = pset1<Packet>(-Scalar(scalar_normalization_offset)); // -24
81 exponent_offset = pselect(is_denormal, padd(exponent_offset, normalization_offset), exponent_offset);
82
83 // Determine exponent and mantissa from normalized_a.
84 exponent = pfrexp_generic_get_biased_exponent(normalized_a);
85 // Zero, Inf and NaN return 'a' unmodified, exponent is zero
86 // (technically the exponent is unspecified for inf/NaN, but GCC/Clang set it to zero)
87 const Scalar scalar_non_finite_exponent = Scalar((ScalarUI(1) << ExponentBits) - ScalarUI(1)); // 255
88 const Packet non_finite_exponent = pset1<Packet>(scalar_non_finite_exponent);
89 const Packet is_zero_or_not_finite = por(pcmp_eq(a, zero), pcmp_eq(exponent, non_finite_exponent));
90 const Packet m = pselect(is_zero_or_not_finite, a, por(pand(normalized_a, sign_mantissa_mask), half));
91 exponent = pselect(is_zero_or_not_finite, zero, padd(exponent, exponent_offset));
92 return m;
93}
94
95// Safely applies ldexp, correctly handles overflows, underflows and denormals.
96// Assumes IEEE floating point format.
97template <typename Packet>
98EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC Packet pldexp_generic(const Packet& a, const Packet& exponent) {
99 // We want to return a * 2^exponent, allowing for all possible integer
100 // exponents without overflowing or underflowing in intermediate
101 // computations.
102 //
103 // Since 'a' and the output can be denormal, the maximum range of 'exponent'
104 // to consider for a float is:
105 // -255-23 -> 255+23
106 // Below -278 any finite float 'a' will become zero, and above +278 any
107 // finite float will become inf, including when 'a' is the smallest possible
108 // denormal.
109 //
110 // Unfortunately, 2^(278) cannot be represented using either one or two
111 // finite normal floats, so we must split the scale factor into at least
112 // three parts. It turns out to be faster to split 'exponent' into four
113 // factors, since [exponent>>2] is much faster to compute that [exponent/3].
114 //
115 // Set e = min(max(exponent, -278), 278);
116 // b = floor(e/4);
117 // out = ((((a * 2^(b)) * 2^(b)) * 2^(b)) * 2^(e-3*b))
118 //
119 // This will avoid any intermediate overflows and correctly handle 0, inf,
120 // NaN cases.
121 typedef typename unpacket_traits<Packet>::integer_packet PacketI;
122 typedef typename unpacket_traits<Packet>::type Scalar;
123 typedef typename unpacket_traits<PacketI>::type ScalarI;
124 static constexpr int TotalBits = sizeof(Scalar) * CHAR_BIT, MantissaBits = numext::numeric_limits<Scalar>::digits - 1,
125 ExponentBits = TotalBits - MantissaBits - 1;
126
127 const Packet max_exponent = pset1<Packet>(Scalar((ScalarI(1) << ExponentBits) + ScalarI(MantissaBits - 1))); // 278
128 const PacketI bias = pset1<PacketI>((ScalarI(1) << (ExponentBits - 1)) - ScalarI(1)); // 127
129 const PacketI e = pcast<Packet, PacketI>(pmin(pmax(exponent, pnegate(max_exponent)), max_exponent));
130 PacketI b = parithmetic_shift_right<2>(e); // floor(e/4);
131 Packet c = preinterpret<Packet>(plogical_shift_left<MantissaBits>(padd(b, bias))); // 2^b
132 Packet out = pmul(pmul(pmul(a, c), c), c); // a * 2^(3b)
133 b = pnmadd(pset1<PacketI>(3), b, e); // e - 3b
134 c = preinterpret<Packet>(plogical_shift_left<MantissaBits>(padd(b, bias))); // 2^(e-3*b)
135 out = pmul(out, c);
136 return out;
137}
138
139// Explicitly multiplies
140// a * (2^e)
141// clamping e to the range
142// [NumTraits<Scalar>::min_exponent()-2, NumTraits<Scalar>::max_exponent()]
143//
144// This is approx 7x faster than pldexp_impl, but will prematurely over/underflow
145// if 2^e doesn't fit into a normal floating-point Scalar.
146//
147// Assumes IEEE floating point format
148template <typename Packet>
149struct pldexp_fast_impl {
150 typedef typename unpacket_traits<Packet>::integer_packet PacketI;
151 typedef typename unpacket_traits<Packet>::type Scalar;
152 typedef typename unpacket_traits<PacketI>::type ScalarI;
153 static constexpr int TotalBits = sizeof(Scalar) * CHAR_BIT, MantissaBits = numext::numeric_limits<Scalar>::digits - 1,
154 ExponentBits = TotalBits - MantissaBits - 1;
155
156 static EIGEN_STRONG_INLINE EIGEN_DEVICE_FUNC Packet run(const Packet& a, const Packet& exponent) {
157 const Packet bias = pset1<Packet>(Scalar((ScalarI(1) << (ExponentBits - 1)) - ScalarI(1))); // 127
158 const Packet limit = pset1<Packet>(Scalar((ScalarI(1) << ExponentBits) - ScalarI(1))); // 255
159 // restrict biased exponent between 0 and 255 for float.
160 const PacketI e = pcast<Packet, PacketI>(pmin(pmax(padd(exponent, bias), pzero(limit)), limit)); // exponent + 127
161 // return a * (2^e)
162 return pmul(a, preinterpret<Packet>(plogical_shift_left<MantissaBits>(e)));
163 }
164};
165
166// Natural or base 2 logarithm.
167// Computes log(x) as log(2^e * m) = C*e + log(m), where the constant C =log(2)
168// and m is in the range [sqrt(1/2),sqrt(2)). In this range, the logarithm can
169// be easily approximated by a polynomial centered on m=1 for stability.
170// TODO(gonnet): Further reduce the interval allowing for lower-degree
171// polynomial interpolants -> ... -> profit!
172template <typename Packet, bool base2>
173EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet plog_impl_float(const Packet _x) {
174 const Packet cst_1 = pset1<Packet>(1.0f);
175 const Packet cst_minus_inf = pset1frombits<Packet>(static_cast<Eigen::numext::uint32_t>(0xff800000u));
176 const Packet cst_pos_inf = pset1frombits<Packet>(static_cast<Eigen::numext::uint32_t>(0x7f800000u));
177
178 const Packet cst_cephes_SQRTHF = pset1<Packet>(0.707106781186547524f);
179 Packet e, x;
180 // extract significant in the range [0.5,1) and exponent
181 x = pfrexp(_x, e);
182
183 // part2: Shift the inputs from the range [0.5,1) to [sqrt(1/2),sqrt(2))
184 // and shift by -1. The values are then centered around 0, which improves
185 // the stability of the polynomial evaluation.
186 // if( x < SQRTHF ) {
187 // e -= 1;
188 // x = x + x - 1.0;
189 // } else { x = x - 1.0; }
190 Packet mask = pcmp_lt(x, cst_cephes_SQRTHF);
191 Packet tmp = pand(x, mask);
192 x = psub(x, cst_1);
193 e = psub(e, pand(cst_1, mask));
194 x = padd(x, tmp);
195
196 // Polynomial coefficients for rational (3,3) r(x) = p(x)/q(x)
197 // approximating log(1+x) on [sqrt(0.5)-1;sqrt(2)-1].
198 const Packet cst_p1 = pset1<Packet>(1.0000000190281136f);
199 const Packet cst_p2 = pset1<Packet>(1.0000000190281063f);
200 const Packet cst_p3 = pset1<Packet>(0.18256296349849254f);
201 const Packet cst_q1 = pset1<Packet>(1.4999999999999927f);
202 const Packet cst_q2 = pset1<Packet>(0.59923249590823520f);
203 const Packet cst_q3 = pset1<Packet>(0.049616247954120038f);
204
205 Packet p = pmadd(x, cst_p3, cst_p2);
206 p = pmadd(x, p, cst_p1);
207 p = pmul(x, p);
208 Packet q = pmadd(x, cst_q3, cst_q2);
209 q = pmadd(x, q, cst_q1);
210 q = pmadd(x, q, cst_1);
211 x = pdiv(p, q);
212
213 // Add the logarithm of the exponent back to the result of the interpolation.
214 if (base2) {
215 const Packet cst_log2e = pset1<Packet>(static_cast<float>(EIGEN_LOG2E));
216 x = pmadd(x, cst_log2e, e);
217 } else {
218 const Packet cst_ln2 = pset1<Packet>(static_cast<float>(EIGEN_LN2));
219 x = pmadd(e, cst_ln2, x);
220 }
221
222 Packet invalid_mask = pcmp_lt_or_nan(_x, pzero(_x));
223 Packet iszero_mask = pcmp_eq(_x, pzero(_x));
224 Packet pos_inf_mask = pcmp_eq(_x, cst_pos_inf);
225 // Filter out invalid inputs, i.e.:
226 // - negative arg will be NAN
227 // - 0 will be -INF
228 // - +INF will be +INF
229 return pselect(iszero_mask, cst_minus_inf, por(pselect(pos_inf_mask, cst_pos_inf, x), invalid_mask));
230}
231
232template <typename Packet>
233EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet plog_float(const Packet _x) {
234 return plog_impl_float<Packet, /* base2 */ false>(_x);
235}
236
237template <typename Packet>
238EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet plog2_float(const Packet _x) {
239 return plog_impl_float<Packet, /* base2 */ true>(_x);
240}
241
242/* Returns the base e (2.718...) or base 2 logarithm of x.
243 * The argument is separated into its exponent and fractional parts.
244 * The logarithm of the fraction in the interval [sqrt(1/2), sqrt(2)],
245 * is approximated by
246 *
247 * log(1+x) = x - 0.5 x**2 + x**3 P(x)/Q(x).
248 *
249 * for more detail see: http://www.netlib.org/cephes/
250 */
251template <typename Packet, bool base2>
252EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet plog_impl_double(const Packet _x) {
253 Packet x = _x;
254
255 const Packet cst_1 = pset1<Packet>(1.0);
256 const Packet cst_neg_half = pset1<Packet>(-0.5);
257 const Packet cst_minus_inf = pset1frombits<Packet>(static_cast<uint64_t>(0xfff0000000000000ull));
258 const Packet cst_pos_inf = pset1frombits<Packet>(static_cast<uint64_t>(0x7ff0000000000000ull));
259
260 // Polynomial Coefficients for log(1+x) = x - x**2/2 + x**3 P(x)/Q(x)
261 // 1/sqrt(2) <= x < sqrt(2)
262 const Packet cst_cephes_SQRTHF = pset1<Packet>(0.70710678118654752440E0);
263 const Packet cst_cephes_log_p0 = pset1<Packet>(1.01875663804580931796E-4);
264 const Packet cst_cephes_log_p1 = pset1<Packet>(4.97494994976747001425E-1);
265 const Packet cst_cephes_log_p2 = pset1<Packet>(4.70579119878881725854E0);
266 const Packet cst_cephes_log_p3 = pset1<Packet>(1.44989225341610930846E1);
267 const Packet cst_cephes_log_p4 = pset1<Packet>(1.79368678507819816313E1);
268 const Packet cst_cephes_log_p5 = pset1<Packet>(7.70838733755885391666E0);
269
270 const Packet cst_cephes_log_q0 = pset1<Packet>(1.0);
271 const Packet cst_cephes_log_q1 = pset1<Packet>(1.12873587189167450590E1);
272 const Packet cst_cephes_log_q2 = pset1<Packet>(4.52279145837532221105E1);
273 const Packet cst_cephes_log_q3 = pset1<Packet>(8.29875266912776603211E1);
274 const Packet cst_cephes_log_q4 = pset1<Packet>(7.11544750618563894466E1);
275 const Packet cst_cephes_log_q5 = pset1<Packet>(2.31251620126765340583E1);
276
277 Packet e;
278 // extract significant in the range [0.5,1) and exponent
279 x = pfrexp(x, e);
280
281 // Shift the inputs from the range [0.5,1) to [sqrt(1/2),sqrt(2))
282 // and shift by -1. The values are then centered around 0, which improves
283 // the stability of the polynomial evaluation.
284 // if( x < SQRTHF ) {
285 // e -= 1;
286 // x = x + x - 1.0;
287 // } else { x = x - 1.0; }
288 Packet mask = pcmp_lt(x, cst_cephes_SQRTHF);
289 Packet tmp = pand(x, mask);
290 x = psub(x, cst_1);
291 e = psub(e, pand(cst_1, mask));
292 x = padd(x, tmp);
293
294 Packet x2 = pmul(x, x);
295 Packet x3 = pmul(x2, x);
296
297 // Evaluate the polynomial approximant , probably to improve instruction-level parallelism.
298 // y = x - 0.5*x^2 + x^3 * polevl( x, P, 5 ) / p1evl( x, Q, 5 ) );
299 Packet y, y1, y_;
300 y = pmadd(cst_cephes_log_p0, x, cst_cephes_log_p1);
301 y1 = pmadd(cst_cephes_log_p3, x, cst_cephes_log_p4);
302 y = pmadd(y, x, cst_cephes_log_p2);
303 y1 = pmadd(y1, x, cst_cephes_log_p5);
304 y_ = pmadd(y, x3, y1);
305
306 y = pmadd(cst_cephes_log_q0, x, cst_cephes_log_q1);
307 y1 = pmadd(cst_cephes_log_q3, x, cst_cephes_log_q4);
308 y = pmadd(y, x, cst_cephes_log_q2);
309 y1 = pmadd(y1, x, cst_cephes_log_q5);
310 y = pmadd(y, x3, y1);
311
312 y_ = pmul(y_, x3);
313 y = pdiv(y_, y);
314
315 y = pmadd(cst_neg_half, x2, y);
316 x = padd(x, y);
317
318 // Add the logarithm of the exponent back to the result of the interpolation.
319 if (base2) {
320 const Packet cst_log2e = pset1<Packet>(static_cast<double>(EIGEN_LOG2E));
321 x = pmadd(x, cst_log2e, e);
322 } else {
323 const Packet cst_ln2 = pset1<Packet>(static_cast<double>(EIGEN_LN2));
324 x = pmadd(e, cst_ln2, x);
325 }
326
327 Packet invalid_mask = pcmp_lt_or_nan(_x, pzero(_x));
328 Packet iszero_mask = pcmp_eq(_x, pzero(_x));
329 Packet pos_inf_mask = pcmp_eq(_x, cst_pos_inf);
330 // Filter out invalid inputs, i.e.:
331 // - negative arg will be NAN
332 // - 0 will be -INF
333 // - +INF will be +INF
334 return pselect(iszero_mask, cst_minus_inf, por(pselect(pos_inf_mask, cst_pos_inf, x), invalid_mask));
335}
336
337template <typename Packet>
338EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet plog_double(const Packet _x) {
339 return plog_impl_double<Packet, /* base2 */ false>(_x);
340}
341
342template <typename Packet>
343EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet plog2_double(const Packet _x) {
344 return plog_impl_double<Packet, /* base2 */ true>(_x);
345}
346
350template <typename Packet>
351Packet generic_plog1p(const Packet& x) {
352 typedef typename unpacket_traits<Packet>::type ScalarType;
353 const Packet one = pset1<Packet>(ScalarType(1));
354 Packet xp1 = padd(x, one);
355 Packet small_mask = pcmp_eq(xp1, one);
356 Packet log1 = plog(xp1);
357 Packet inf_mask = pcmp_eq(xp1, log1);
358 Packet log_large = pmul(x, pdiv(log1, psub(xp1, one)));
359 return pselect(por(small_mask, inf_mask), x, log_large);
360}
361
365template <typename Packet>
366Packet generic_expm1(const Packet& x) {
367 typedef typename unpacket_traits<Packet>::type ScalarType;
368 const Packet one = pset1<Packet>(ScalarType(1));
369 const Packet neg_one = pset1<Packet>(ScalarType(-1));
370 Packet u = pexp(x);
371 Packet one_mask = pcmp_eq(u, one);
372 Packet u_minus_one = psub(u, one);
373 Packet neg_one_mask = pcmp_eq(u_minus_one, neg_one);
374 Packet logu = plog(u);
375 // The following comparison is to catch the case where
376 // exp(x) = +inf. It is written in this way to avoid having
377 // to form the constant +inf, which depends on the packet
378 // type.
379 Packet pos_inf_mask = pcmp_eq(logu, u);
380 Packet expm1 = pmul(u_minus_one, pdiv(x, logu));
381 expm1 = pselect(pos_inf_mask, u, expm1);
382 return pselect(one_mask, x, pselect(neg_one_mask, neg_one, expm1));
383}
384
385// Exponential function. Works by writing "x = m*log(2) + r" where
386// "m = floor(x/log(2)+1/2)" and "r" is the remainder. The result is then
387// "exp(x) = 2^m*exp(r)" where exp(r) is in the range [-1,1).
388// exp(r) is computed using a 6th order minimax polynomial approximation.
389template <typename Packet>
390EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet pexp_float(const Packet _x) {
391 const Packet cst_zero = pset1<Packet>(0.0f);
392 const Packet cst_one = pset1<Packet>(1.0f);
393 const Packet cst_half = pset1<Packet>(0.5f);
394 const Packet cst_exp_hi = pset1<Packet>(88.723f);
395 const Packet cst_exp_lo = pset1<Packet>(-104.f);
396
397 const Packet cst_cephes_LOG2EF = pset1<Packet>(1.44269504088896341f);
398 const Packet cst_p2 = pset1<Packet>(0.49999988079071044921875f);
399 const Packet cst_p3 = pset1<Packet>(0.16666518151760101318359375f);
400 const Packet cst_p4 = pset1<Packet>(4.166965186595916748046875e-2f);
401 const Packet cst_p5 = pset1<Packet>(8.36894474923610687255859375e-3f);
402 const Packet cst_p6 = pset1<Packet>(1.37449637986719608306884765625e-3f);
403
404 // Clamp x.
405 Packet zero_mask = pcmp_lt(_x, cst_exp_lo);
406 Packet x = pmin(_x, cst_exp_hi);
407
408 // Express exp(x) as exp(m*ln(2) + r), start by extracting
409 // m = floor(x/ln(2) + 0.5).
410 Packet m = pfloor(pmadd(x, cst_cephes_LOG2EF, cst_half));
411
412 // Get r = x - m*ln(2). If no FMA instructions are available, m*ln(2) is
413 // subtracted out in two parts, m*C1+m*C2 = m*ln(2), to avoid accumulating
414 // truncation errors.
415 const Packet cst_cephes_exp_C1 = pset1<Packet>(-0.693359375f);
416 const Packet cst_cephes_exp_C2 = pset1<Packet>(2.12194440e-4f);
417 Packet r = pmadd(m, cst_cephes_exp_C1, x);
418 r = pmadd(m, cst_cephes_exp_C2, r);
419
420 // Evaluate the 6th order polynomial approximation to exp(r)
421 // with r in the interval [-ln(2)/2;ln(2)/2].
422 const Packet r2 = pmul(r, r);
423 Packet p_even = pmadd(r2, cst_p6, cst_p4);
424 const Packet p_odd = pmadd(r2, cst_p5, cst_p3);
425 p_even = pmadd(r2, p_even, cst_p2);
426 const Packet p_low = padd(r, cst_one);
427 Packet y = pmadd(r, p_odd, p_even);
428 y = pmadd(r2, y, p_low);
429
430 // Return 2^m * exp(r).
431 // TODO: replace pldexp with faster implementation since y in [-1, 1).
432 return pselect(zero_mask, cst_zero, pmax(pldexp(y, m), _x));
433}
434
435template <typename Packet>
436EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet pexp_double(const Packet _x) {
437 Packet x = _x;
438 const Packet cst_zero = pset1<Packet>(0.0);
439 const Packet cst_1 = pset1<Packet>(1.0);
440 const Packet cst_2 = pset1<Packet>(2.0);
441 const Packet cst_half = pset1<Packet>(0.5);
442
443 const Packet cst_exp_hi = pset1<Packet>(709.784);
444 const Packet cst_exp_lo = pset1<Packet>(-709.784);
445
446 const Packet cst_cephes_LOG2EF = pset1<Packet>(1.4426950408889634073599);
447 const Packet cst_cephes_exp_p0 = pset1<Packet>(1.26177193074810590878e-4);
448 const Packet cst_cephes_exp_p1 = pset1<Packet>(3.02994407707441961300e-2);
449 const Packet cst_cephes_exp_p2 = pset1<Packet>(9.99999999999999999910e-1);
450 const Packet cst_cephes_exp_q0 = pset1<Packet>(3.00198505138664455042e-6);
451 const Packet cst_cephes_exp_q1 = pset1<Packet>(2.52448340349684104192e-3);
452 const Packet cst_cephes_exp_q2 = pset1<Packet>(2.27265548208155028766e-1);
453 const Packet cst_cephes_exp_q3 = pset1<Packet>(2.00000000000000000009e0);
454 const Packet cst_cephes_exp_C1 = pset1<Packet>(0.693145751953125);
455 const Packet cst_cephes_exp_C2 = pset1<Packet>(1.42860682030941723212e-6);
456
457 Packet tmp, fx;
458
459 // clamp x
460 Packet zero_mask = pcmp_lt(_x, cst_exp_lo);
461 x = pmin(x, cst_exp_hi);
462 // Express exp(x) as exp(g + n*log(2)).
463 fx = pmadd(cst_cephes_LOG2EF, x, cst_half);
464
465 // Get the integer modulus of log(2), i.e. the "n" described above.
466 fx = pfloor(fx);
467
468 // Get the remainder modulo log(2), i.e. the "g" described above. Subtract
469 // n*log(2) out in two steps, i.e. n*C1 + n*C2, C1+C2=log2 to get the last
470 // digits right.
471 tmp = pmul(fx, cst_cephes_exp_C1);
472 Packet z = pmul(fx, cst_cephes_exp_C2);
473 x = psub(x, tmp);
474 x = psub(x, z);
475
476 Packet x2 = pmul(x, x);
477
478 // Evaluate the numerator polynomial of the rational interpolant.
479 Packet px = cst_cephes_exp_p0;
480 px = pmadd(px, x2, cst_cephes_exp_p1);
481 px = pmadd(px, x2, cst_cephes_exp_p2);
482 px = pmul(px, x);
483
484 // Evaluate the denominator polynomial of the rational interpolant.
485 Packet qx = cst_cephes_exp_q0;
486 qx = pmadd(qx, x2, cst_cephes_exp_q1);
487 qx = pmadd(qx, x2, cst_cephes_exp_q2);
488 qx = pmadd(qx, x2, cst_cephes_exp_q3);
489
490 // I don't really get this bit, copied from the SSE2 routines, so...
491 // TODO(gonnet): Figure out what is going on here, perhaps find a better
492 // rational interpolant?
493 x = pdiv(px, psub(qx, px));
494 x = pmadd(cst_2, x, cst_1);
495
496 // Construct the result 2^n * exp(g) = e * x. The max is used to catch
497 // non-finite values in the input.
498 // TODO: replace pldexp with faster implementation since x in [-1, 1).
499 return pselect(zero_mask, cst_zero, pmax(pldexp(x, fx), _x));
500}
501
502// The following code is inspired by the following stack-overflow answer:
503// https://stackoverflow.com/questions/30463616/payne-hanek-algorithm-implementation-in-c/30465751#30465751
504// It has been largely optimized:
505// - By-pass calls to frexp.
506// - Aligned loads of required 96 bits of 2/pi. This is accomplished by
507// (1) balancing the mantissa and exponent to the required bits of 2/pi are
508// aligned on 8-bits, and (2) replicating the storage of the bits of 2/pi.
509// - Avoid a branch in rounding and extraction of the remaining fractional part.
510// Overall, I measured a speed up higher than x2 on x86-64.
511inline float trig_reduce_huge(float xf, Eigen::numext::int32_t* quadrant) {
512 using Eigen::numext::int32_t;
513 using Eigen::numext::int64_t;
514 using Eigen::numext::uint32_t;
515 using Eigen::numext::uint64_t;
516
517 const double pio2_62 = 3.4061215800865545e-19; // pi/2 * 2^-62
518 const uint64_t zero_dot_five = uint64_t(1) << 61; // 0.5 in 2.62-bit fixed-point format
519
520 // 192 bits of 2/pi for Payne-Hanek reduction
521 // Bits are introduced by packet of 8 to enable aligned reads.
522 static const uint32_t two_over_pi[] = {
523 0x00000028, 0x000028be, 0x0028be60, 0x28be60db, 0xbe60db93, 0x60db9391, 0xdb939105, 0x9391054a, 0x91054a7f,
524 0x054a7f09, 0x4a7f09d5, 0x7f09d5f4, 0x09d5f47d, 0xd5f47d4d, 0xf47d4d37, 0x7d4d3770, 0x4d377036, 0x377036d8,
525 0x7036d8a5, 0x36d8a566, 0xd8a5664f, 0xa5664f10, 0x664f10e4, 0x4f10e410, 0x10e41000, 0xe4100000};
526
527 uint32_t xi = numext::bit_cast<uint32_t>(xf);
528 // Below, -118 = -126 + 8.
529 // -126 is to get the exponent,
530 // +8 is to enable alignment of 2/pi's bits on 8 bits.
531 // This is possible because the fractional part of x as only 24 meaningful bits.
532 uint32_t e = (xi >> 23) - 118;
533 // Extract the mantissa and shift it to align it wrt the exponent
534 xi = ((xi & 0x007fffffu) | 0x00800000u) << (e & 0x7);
535
536 uint32_t i = e >> 3;
537 uint32_t twoopi_1 = two_over_pi[i - 1];
538 uint32_t twoopi_2 = two_over_pi[i + 3];
539 uint32_t twoopi_3 = two_over_pi[i + 7];
540
541 // Compute x * 2/pi in 2.62-bit fixed-point format.
542 uint64_t p;
543 p = uint64_t(xi) * twoopi_3;
544 p = uint64_t(xi) * twoopi_2 + (p >> 32);
545 p = (uint64_t(xi * twoopi_1) << 32) + p;
546
547 // Round to nearest: add 0.5 and extract integral part.
548 uint64_t q = (p + zero_dot_five) >> 62;
549 *quadrant = int(q);
550 // Now it remains to compute "r = x - q*pi/2" with high accuracy,
551 // since we have p=x/(pi/2) with high accuracy, we can more efficiently compute r as:
552 // r = (p-q)*pi/2,
553 // where the product can be be carried out with sufficient accuracy using double precision.
554 p -= q << 62;
555 return float(double(int64_t(p)) * pio2_62);
556}
557
558template <bool ComputeSine, typename Packet, bool ComputeBoth = false>
559EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS
560#if EIGEN_COMP_GNUC_STRICT
561 __attribute__((optimize("-fno-unsafe-math-optimizations")))
562#endif
563 Packet
564 psincos_float(const Packet& _x) {
565 typedef typename unpacket_traits<Packet>::integer_packet PacketI;
566
567 const Packet cst_2oPI = pset1<Packet>(0.636619746685028076171875f); // 2/PI
568 const Packet cst_rounding_magic = pset1<Packet>(12582912); // 2^23 for rounding
569 const PacketI csti_1 = pset1<PacketI>(1);
570 const Packet cst_sign_mask = pset1frombits<Packet>(static_cast<Eigen::numext::uint32_t>(0x80000000u));
571
572 Packet x = pabs(_x);
573
574 // Scale x by 2/Pi to find x's octant.
575 Packet y = pmul(x, cst_2oPI);
576
577 // Rounding trick to find nearest integer:
578 Packet y_round = padd(y, cst_rounding_magic);
579 EIGEN_OPTIMIZATION_BARRIER(y_round)
580 PacketI y_int = preinterpret<PacketI>(y_round); // last 23 digits represent integer (if abs(x)<2^24)
581 y = psub(y_round, cst_rounding_magic); // nearest integer to x * (2/pi)
582
583// Subtract y * Pi/2 to reduce x to the interval -Pi/4 <= x <= +Pi/4
584// using "Extended precision modular arithmetic"
585#if defined(EIGEN_VECTORIZE_FMA)
586 // This version requires true FMA for high accuracy.
587 // It provides a max error of 1ULP up to (with absolute_error < 5.9605e-08):
588 const float huge_th = ComputeSine ? 117435.992f : 71476.0625f;
589 x = pmadd(y, pset1<Packet>(-1.57079601287841796875f), x);
590 x = pmadd(y, pset1<Packet>(-3.1391647326017846353352069854736328125e-07f), x);
591 x = pmadd(y, pset1<Packet>(-5.390302529957764765544681040410068817436695098876953125e-15f), x);
592#else
593 // Without true FMA, the previous set of coefficients maintain 1ULP accuracy
594 // up to x<15.7 (for sin), but accuracy is immediately lost for x>15.7.
595 // We thus use one more iteration to maintain 2ULPs up to reasonably large inputs.
596
597 // The following set of coefficients maintain 1ULP up to 9.43 and 14.16 for sin and cos respectively.
598 // and 2 ULP up to:
599 const float huge_th = ComputeSine ? 25966.f : 18838.f;
600 x = pmadd(y, pset1<Packet>(-1.5703125), x); // = 0xbfc90000
601 EIGEN_OPTIMIZATION_BARRIER(x)
602 x = pmadd(y, pset1<Packet>(-0.000483989715576171875), x); // = 0xb9fdc000
603 EIGEN_OPTIMIZATION_BARRIER(x)
604 x = pmadd(y, pset1<Packet>(1.62865035235881805419921875e-07), x); // = 0x342ee000
605 x = pmadd(y, pset1<Packet>(5.5644315544167710640977020375430583953857421875e-11), x); // = 0x2e74b9ee
606
607// For the record, the following set of coefficients maintain 2ULP up
608// to a slightly larger range:
609// const float huge_th = ComputeSine ? 51981.f : 39086.125f;
610// but it slightly fails to maintain 1ULP for two values of sin below pi.
611// x = pmadd(y, pset1<Packet>(-3.140625/2.), x);
612// x = pmadd(y, pset1<Packet>(-0.00048351287841796875), x);
613// x = pmadd(y, pset1<Packet>(-3.13855707645416259765625e-07), x);
614// x = pmadd(y, pset1<Packet>(-6.0771006282767103812147979624569416046142578125e-11), x);
615
616// For the record, with only 3 iterations it is possible to maintain
617// 1 ULP up to 3PI (maybe more) and 2ULP up to 255.
618// The coefficients are: 0xbfc90f80, 0xb7354480, 0x2e74b9ee
619#endif
620
621 if (predux_any(pcmp_le(pset1<Packet>(huge_th), pabs(_x)))) {
622 const int PacketSize = unpacket_traits<Packet>::size;
623 EIGEN_ALIGN_TO_BOUNDARY(sizeof(Packet)) float vals[PacketSize];
624 EIGEN_ALIGN_TO_BOUNDARY(sizeof(Packet)) float x_cpy[PacketSize];
625 EIGEN_ALIGN_TO_BOUNDARY(sizeof(Packet)) Eigen::numext::int32_t y_int2[PacketSize];
626 pstoreu(vals, pabs(_x));
627 pstoreu(x_cpy, x);
628 pstoreu(y_int2, y_int);
629 for (int k = 0; k < PacketSize; ++k) {
630 float val = vals[k];
631 if (val >= huge_th && (numext::isfinite)(val)) x_cpy[k] = trig_reduce_huge(val, &y_int2[k]);
632 }
633 x = ploadu<Packet>(x_cpy);
634 y_int = ploadu<PacketI>(y_int2);
635 }
636
637 // Compute the sign to apply to the polynomial.
638 // sin: sign = second_bit(y_int) xor signbit(_x)
639 // cos: sign = second_bit(y_int+1)
640 Packet sign_bit = ComputeSine ? pxor(_x, preinterpret<Packet>(plogical_shift_left<30>(y_int)))
641 : preinterpret<Packet>(plogical_shift_left<30>(padd(y_int, csti_1)));
642 sign_bit = pand(sign_bit, cst_sign_mask); // clear all but left most bit
643
644 // Get the polynomial selection mask from the second bit of y_int
645 // We'll calculate both (sin and cos) polynomials and then select from the two.
646 Packet poly_mask = preinterpret<Packet>(pcmp_eq(pand(y_int, csti_1), pzero(y_int)));
647
648 Packet x2 = pmul(x, x);
649
650 // Evaluate the cos(x) polynomial. (-Pi/4 <= x <= Pi/4)
651 Packet y1 = pset1<Packet>(2.4372266125283204019069671630859375e-05f);
652 y1 = pmadd(y1, x2, pset1<Packet>(-0.00138865201734006404876708984375f));
653 y1 = pmadd(y1, x2, pset1<Packet>(0.041666619479656219482421875f));
654 y1 = pmadd(y1, x2, pset1<Packet>(-0.5f));
655 y1 = pmadd(y1, x2, pset1<Packet>(1.f));
656
657 // Evaluate the sin(x) polynomial. (Pi/4 <= x <= Pi/4)
658 // octave/matlab code to compute those coefficients:
659 // x = (0:0.0001:pi/4)';
660 // A = [x.^3 x.^5 x.^7];
661 // w = ((1.-(x/(pi/4)).^2).^5)*2000+1; # weights trading relative accuracy
662 // c = (A'*diag(w)*A)\‍(A'*diag(w)*(sin(x)-x)); # weighted LS, linear coeff forced to 1
663 // printf('%.64f\n %.64f\n%.64f\n', c(3), c(2), c(1))
664 //
665 Packet y2 = pset1<Packet>(-0.0001959234114083702898469196984621021329076029360294342041015625f);
666 y2 = pmadd(y2, x2, pset1<Packet>(0.0083326873655616851693794799871284340042620897293090820312500000f));
667 y2 = pmadd(y2, x2, pset1<Packet>(-0.1666666203982298255503735617821803316473960876464843750000000000f));
668 y2 = pmul(y2, x2);
669 y2 = pmadd(y2, x, x);
670
671 // Select the correct result from the two polynomials.
672 if (ComputeBoth) {
673 Packet peven = peven_mask(x);
674 Packet ysin = pselect(poly_mask, y2, y1);
675 Packet ycos = pselect(poly_mask, y1, y2);
676 Packet sign_bit_sin = pxor(_x, preinterpret<Packet>(plogical_shift_left<30>(y_int)));
677 Packet sign_bit_cos = preinterpret<Packet>(plogical_shift_left<30>(padd(y_int, csti_1)));
678 sign_bit_sin = pand(sign_bit_sin, cst_sign_mask); // clear all but left most bit
679 sign_bit_cos = pand(sign_bit_cos, cst_sign_mask); // clear all but left most bit
680 y = pselect(peven, pxor(ysin, sign_bit_sin), pxor(ycos, sign_bit_cos));
681 } else {
682 y = ComputeSine ? pselect(poly_mask, y2, y1) : pselect(poly_mask, y1, y2);
683 y = pxor(y, sign_bit);
684 }
685 // Update the sign and filter huge inputs
686 return y;
687}
688
689template <typename Packet>
690EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet psin_float(const Packet& x) {
691 return psincos_float<true>(x);
692}
693
694template <typename Packet>
695EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet pcos_float(const Packet& x) {
696 return psincos_float<false>(x);
697}
698
699// Trigonometric argument reduction for double for inputs smaller than 15.
700// Reduces trigonometric arguments for double inputs where x < 15. Given an argument x and its corresponding quadrant
701// count n, the function computes and returns the reduced argument t such that x = n * pi/2 + t.
702template <typename Packet>
703Packet trig_reduce_small_double(const Packet& x, const Packet& q) {
704 // Pi/2 split into 2 values
705 const Packet cst_pio2_a = pset1<Packet>(-1.570796325802803);
706 const Packet cst_pio2_b = pset1<Packet>(-9.920935184482005e-10);
707
708 Packet t;
709 t = pmadd(cst_pio2_a, q, x);
710 t = pmadd(cst_pio2_b, q, t);
711 return t;
712}
713
714// Trigonometric argument reduction for double for inputs smaller than 1e14.
715// Reduces trigonometric arguments for double inputs where x < 1e14. Given an argument x and its corresponding quadrant
716// count n, the function computes and returns the reduced argument t such that x = n * pi/2 + t.
717template <typename Packet>
718Packet trig_reduce_medium_double(const Packet& x, const Packet& q_high, const Packet& q_low) {
719 // Pi/2 split into 4 values
720 const Packet cst_pio2_a = pset1<Packet>(-1.570796325802803);
721 const Packet cst_pio2_b = pset1<Packet>(-9.920935184482005e-10);
722 const Packet cst_pio2_c = pset1<Packet>(-6.123234014771656e-17);
723 const Packet cst_pio2_d = pset1<Packet>(1.903488962019325e-25);
724
725 Packet t;
726 t = pmadd(cst_pio2_a, q_high, x);
727 t = pmadd(cst_pio2_a, q_low, t);
728 t = pmadd(cst_pio2_b, q_high, t);
729 t = pmadd(cst_pio2_b, q_low, t);
730 t = pmadd(cst_pio2_c, q_high, t);
731 t = pmadd(cst_pio2_c, q_low, t);
732 t = pmadd(cst_pio2_d, padd(q_low, q_high), t);
733 return t;
734}
735
736template <bool ComputeSine, typename Packet, bool ComputeBoth = false>
737EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS
738#if EIGEN_COMP_GNUC_STRICT
739 __attribute__((optimize("-fno-unsafe-math-optimizations")))
740#endif
741 Packet
742 psincos_double(const Packet& x) {
743 typedef typename unpacket_traits<Packet>::integer_packet PacketI;
744 typedef typename unpacket_traits<PacketI>::type ScalarI;
745
746 const Packet cst_sign_mask = pset1frombits<Packet>(static_cast<Eigen::numext::uint64_t>(0x8000000000000000u));
747
748 // If the argument is smaller than this value, use a simpler argument reduction
749 const double small_th = 15;
750 // If the argument is bigger than this value, use the non-vectorized std version
751 const double huge_th = 1e14;
752
753 const Packet cst_2oPI = pset1<Packet>(0.63661977236758134307553505349006); // 2/PI
754 // Integer Packet constants
755 const PacketI cst_one = pset1<PacketI>(ScalarI(1));
756 // Constant for splitting
757 const Packet cst_split = pset1<Packet>(1 << 24);
758
759 Packet x_abs = pabs(x);
760
761 // Scale x by 2/Pi
762 PacketI q_int;
763 Packet s;
764
765 // TODO Implement huge angle argument reduction
766 if (EIGEN_PREDICT_FALSE(predux_any(pcmp_le(pset1<Packet>(small_th), x_abs)))) {
767 Packet q_high = pmul(pfloor(pmul(x_abs, pdiv(cst_2oPI, cst_split))), cst_split);
768 Packet q_low_noround = psub(pmul(x_abs, cst_2oPI), q_high);
769 q_int = pcast<Packet, PacketI>(padd(q_low_noround, pset1<Packet>(0.5)));
770 Packet q_low = pcast<PacketI, Packet>(q_int);
771 s = trig_reduce_medium_double(x_abs, q_high, q_low);
772 } else {
773 Packet qval_noround = pmul(x_abs, cst_2oPI);
774 q_int = pcast<Packet, PacketI>(padd(qval_noround, pset1<Packet>(0.5)));
775 Packet q = pcast<PacketI, Packet>(q_int);
776 s = trig_reduce_small_double(x_abs, q);
777 }
778
779 // All the upcoming approximating polynomials have even exponents
780 Packet ss = pmul(s, s);
781
782 // Padé approximant of cos(x)
783 // Assuring < 1 ULP error on the interval [-pi/4, pi/4]
784 // cos(x) ~= (80737373*x^8 - 13853547000*x^6 + 727718024880*x^4 - 11275015752000*x^2 + 23594700729600)/(147173*x^8 +
785 // 39328920*x^6 + 5772800880*x^4 + 522334612800*x^2 + 23594700729600)
786 // MATLAB code to compute those coefficients:
787 // syms x;
788 // cosf = @(x) cos(x);
789 // pade_cosf = pade(cosf(x), x, 0, 'Order', 8)
790 Packet sc1_num = pmadd(ss, pset1<Packet>(80737373), pset1<Packet>(-13853547000));
791 Packet sc2_num = pmadd(sc1_num, ss, pset1<Packet>(727718024880));
792 Packet sc3_num = pmadd(sc2_num, ss, pset1<Packet>(-11275015752000));
793 Packet sc4_num = pmadd(sc3_num, ss, pset1<Packet>(23594700729600));
794 Packet sc1_denum = pmadd(ss, pset1<Packet>(147173), pset1<Packet>(39328920));
795 Packet sc2_denum = pmadd(sc1_denum, ss, pset1<Packet>(5772800880));
796 Packet sc3_denum = pmadd(sc2_denum, ss, pset1<Packet>(522334612800));
797 Packet sc4_denum = pmadd(sc3_denum, ss, pset1<Packet>(23594700729600));
798 Packet scos = pdiv(sc4_num, sc4_denum);
799
800 // Padé approximant of sin(x)
801 // Assuring < 1 ULP error on the interval [-pi/4, pi/4]
802 // sin(x) ~= (x*(4585922449*x^8 - 1066023933480*x^6 + 83284044283440*x^4 - 2303682236856000*x^2 +
803 // 15605159573203200))/(45*(1029037*x^8 + 345207016*x^6 + 61570292784*x^4 + 6603948711360*x^2 + 346781323848960))
804 // MATLAB code to compute those coefficients:
805 // syms x;
806 // sinf = @(x) sin(x);
807 // pade_sinf = pade(sinf(x), x, 0, 'Order', 8, 'OrderMode', 'relative')
808 Packet ss1_num = pmadd(ss, pset1<Packet>(4585922449), pset1<Packet>(-1066023933480));
809 Packet ss2_num = pmadd(ss1_num, ss, pset1<Packet>(83284044283440));
810 Packet ss3_num = pmadd(ss2_num, ss, pset1<Packet>(-2303682236856000));
811 Packet ss4_num = pmadd(ss3_num, ss, pset1<Packet>(15605159573203200));
812 Packet ss1_denum = pmadd(ss, pset1<Packet>(1029037), pset1<Packet>(345207016));
813 Packet ss2_denum = pmadd(ss1_denum, ss, pset1<Packet>(61570292784));
814 Packet ss3_denum = pmadd(ss2_denum, ss, pset1<Packet>(6603948711360));
815 Packet ss4_denum = pmadd(ss3_denum, ss, pset1<Packet>(346781323848960));
816 Packet ssin = pdiv(pmul(s, ss4_num), pmul(pset1<Packet>(45), ss4_denum));
817
818 Packet poly_mask = preinterpret<Packet>(pcmp_eq(pand(q_int, cst_one), pzero(q_int)));
819
820 Packet sign_sin = pxor(x, preinterpret<Packet>(plogical_shift_left<62>(q_int)));
821 Packet sign_cos = preinterpret<Packet>(plogical_shift_left<62>(padd(q_int, cst_one)));
822 Packet sign_bit, sFinalRes;
823 if (ComputeBoth) {
824 Packet peven = peven_mask(x);
825 sign_bit = pselect((s), sign_sin, sign_cos);
826 sFinalRes = pselect(pxor(peven, poly_mask), ssin, scos);
827 } else {
828 sign_bit = ComputeSine ? sign_sin : sign_cos;
829 sFinalRes = ComputeSine ? pselect(poly_mask, ssin, scos) : pselect(poly_mask, scos, ssin);
830 }
831 sign_bit = pand(sign_bit, cst_sign_mask); // clear all but left most bit
832 sFinalRes = pxor(sFinalRes, sign_bit);
833
834 // If the inputs values are higher than that a value that the argument reduction can currently address, compute them
835 // using std::sin and std::cos
836 // TODO Remove it when huge angle argument reduction is implemented
837 if (EIGEN_PREDICT_FALSE(predux_any(pcmp_le(pset1<Packet>(huge_th), x_abs)))) {
838 const int PacketSize = unpacket_traits<Packet>::size;
839 EIGEN_ALIGN_TO_BOUNDARY(sizeof(Packet)) double sincos_vals[PacketSize];
840 EIGEN_ALIGN_TO_BOUNDARY(sizeof(Packet)) double x_cpy[PacketSize];
841 pstoreu(x_cpy, x);
842 pstoreu(sincos_vals, sFinalRes);
843 for (int k = 0; k < PacketSize; ++k) {
844 double val = x_cpy[k];
845 if (std::abs(val) > huge_th && (numext::isfinite)(val)) {
846 if (ComputeBoth)
847 sincos_vals[k] = k % 2 == 0 ? std::sin(val) : std::cos(val);
848 else
849 sincos_vals[k] = ComputeSine ? std::sin(val) : std::cos(val);
850 }
851 }
852 sFinalRes = ploadu<Packet>(sincos_vals);
853 }
854 return sFinalRes;
855}
856
857template <typename Packet>
858EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet psin_double(const Packet& x) {
859 return psincos_double<true>(x);
860}
861
862template <typename Packet>
863EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet pcos_double(const Packet& x) {
864 return psincos_double<false>(x);
865}
866
867// Generic implementation of acos(x).
868template <typename Packet>
869EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet pacos_float(const Packet& x_in) {
870 typedef typename unpacket_traits<Packet>::type Scalar;
871 static_assert(std::is_same<Scalar, float>::value, "Scalar type must be float");
872
873 const Packet cst_one = pset1<Packet>(Scalar(1));
874 const Packet cst_pi = pset1<Packet>(Scalar(EIGEN_PI));
875 const Packet p6 = pset1<Packet>(Scalar(2.36423197202384471893310546875e-3));
876 const Packet p5 = pset1<Packet>(Scalar(-1.1368644423782825469970703125e-2));
877 const Packet p4 = pset1<Packet>(Scalar(2.717843465507030487060546875e-2));
878 const Packet p3 = pset1<Packet>(Scalar(-4.8969544470310211181640625e-2));
879 const Packet p2 = pset1<Packet>(Scalar(8.8804088532924652099609375e-2));
880 const Packet p1 = pset1<Packet>(Scalar(-0.214591205120086669921875));
881 const Packet p0 = pset1<Packet>(Scalar(1.57079637050628662109375));
882
883 // For x in [0:1], we approximate acos(x)/sqrt(1-x), which is a smooth
884 // function, by a 6'th order polynomial.
885 // For x in [-1:0) we use that acos(-x) = pi - acos(x).
886 const Packet neg_mask = psignbit(x_in);
887 const Packet abs_x = pabs(x_in);
888
889 // Evaluate the polynomial using Horner's rule:
890 // P(x) = p0 + x * (p1 + x * (p2 + ... (p5 + x * p6)) ... ) .
891 // We evaluate even and odd terms independently to increase
892 // instruction level parallelism.
893 Packet x2 = pmul(x_in, x_in);
894 Packet p_even = pmadd(p6, x2, p4);
895 Packet p_odd = pmadd(p5, x2, p3);
896 p_even = pmadd(p_even, x2, p2);
897 p_odd = pmadd(p_odd, x2, p1);
898 p_even = pmadd(p_even, x2, p0);
899 Packet p = pmadd(p_odd, abs_x, p_even);
900
901 // The polynomial approximates acos(x)/sqrt(1-x), so
902 // multiply by sqrt(1-x) to get acos(x).
903 // Conveniently returns NaN for arguments outside [-1:1].
904 Packet denom = psqrt(psub(cst_one, abs_x));
905 Packet result = pmul(denom, p);
906 // Undo mapping for negative arguments.
907 return pselect(neg_mask, psub(cst_pi, result), result);
908}
909
910// Generic implementation of asin(x).
911template <typename Packet>
912EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet pasin_float(const Packet& x_in) {
913 typedef typename unpacket_traits<Packet>::type Scalar;
914 static_assert(std::is_same<Scalar, float>::value, "Scalar type must be float");
915
916 constexpr float kPiOverTwo = static_cast<float>(EIGEN_PI / 2);
917
918 const Packet cst_half = pset1<Packet>(0.5f);
919 const Packet cst_one = pset1<Packet>(1.0f);
920 const Packet cst_two = pset1<Packet>(2.0f);
921 const Packet cst_pi_over_two = pset1<Packet>(kPiOverTwo);
922 // For |x| < 0.5 approximate asin(x)/x by an 8th order polynomial with
923 // even terms only.
924 const Packet p9 = pset1<Packet>(5.08838854730129241943359375e-2f);
925 const Packet p7 = pset1<Packet>(3.95139865577220916748046875e-2f);
926 const Packet p5 = pset1<Packet>(7.550220191478729248046875e-2f);
927 const Packet p3 = pset1<Packet>(0.16664917767047882080078125f);
928 const Packet p1 = pset1<Packet>(1.00000011920928955078125f);
929
930 const Packet abs_x = pabs(x_in);
931 const Packet sign_mask = pandnot(x_in, abs_x);
932 const Packet invalid_mask = pcmp_lt(cst_one, abs_x);
933
934 // For arguments |x| > 0.5, we map x back to [0:0.5] using
935 // the transformation x_large = sqrt(0.5*(1-x)), and use the
936 // identity
937 // asin(x) = pi/2 - 2 * asin( sqrt( 0.5 * (1 - x)))
938
939 const Packet x_large = psqrt(pnmadd(cst_half, abs_x, cst_half));
940 const Packet large_mask = pcmp_lt(cst_half, abs_x);
941 const Packet x = pselect(large_mask, x_large, abs_x);
942 const Packet x2 = pmul(x, x);
943
944 // Compute polynomial.
945 // x * (p1 + x^2*(p3 + x^2*(p5 + x^2*(p7 + x^2*p9))))
946
947 Packet p = pmadd(p9, x2, p7);
948 p = pmadd(p, x2, p5);
949 p = pmadd(p, x2, p3);
950 p = pmadd(p, x2, p1);
951 p = pmul(p, x);
952
953 const Packet p_large = pnmadd(cst_two, p, cst_pi_over_two);
954 p = pselect(large_mask, p_large, p);
955 // Flip the sign for negative arguments.
956 p = pxor(p, sign_mask);
957 // Return NaN for arguments outside [-1:1].
958 return por(invalid_mask, p);
959}
960
961// Computes elementwise atan(x) for x in [-1:1] with 2 ulp accuracy.
962template <typename Packet>
963EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet patan_reduced_float(const Packet& x) {
964 const Packet q0 = pset1<Packet>(-0.3333314359188079833984375f);
965 const Packet q2 = pset1<Packet>(0.19993579387664794921875f);
966 const Packet q4 = pset1<Packet>(-0.14209578931331634521484375f);
967 const Packet q6 = pset1<Packet>(0.1066047251224517822265625f);
968 const Packet q8 = pset1<Packet>(-7.5408883392810821533203125e-2f);
969 const Packet q10 = pset1<Packet>(4.3082617223262786865234375e-2f);
970 const Packet q12 = pset1<Packet>(-1.62907354533672332763671875e-2f);
971 const Packet q14 = pset1<Packet>(2.90188402868807315826416015625e-3f);
972
973 // Approximate atan(x) by a polynomial of the form
974 // P(x) = x + x^3 * Q(x^2),
975 // where Q(x^2) is a 7th order polynomial in x^2.
976 // We evaluate even and odd terms in x^2 in parallel
977 // to take advantage of instruction level parallelism
978 // and hardware with multiple FMA units.
979
980 // note: if x == -0, this returns +0
981 const Packet x2 = pmul(x, x);
982 const Packet x4 = pmul(x2, x2);
983 Packet q_odd = pmadd(q14, x4, q10);
984 Packet q_even = pmadd(q12, x4, q8);
985 q_odd = pmadd(q_odd, x4, q6);
986 q_even = pmadd(q_even, x4, q4);
987 q_odd = pmadd(q_odd, x4, q2);
988 q_even = pmadd(q_even, x4, q0);
989 const Packet q = pmadd(q_odd, x2, q_even);
990 return pmadd(q, pmul(x, x2), x);
991}
992
993template <typename Packet>
994EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet patan_float(const Packet& x_in) {
995 typedef typename unpacket_traits<Packet>::type Scalar;
996 static_assert(std::is_same<Scalar, float>::value, "Scalar type must be float");
997
998 constexpr float kPiOverTwo = static_cast<float>(EIGEN_PI / 2);
999
1000 const Packet cst_signmask = pset1<Packet>(-0.0f);
1001 const Packet cst_one = pset1<Packet>(1.0f);
1002 const Packet cst_pi_over_two = pset1<Packet>(kPiOverTwo);
1003
1004 // "Large": For |x| > 1, use atan(1/x) = sign(x)*pi/2 - atan(x).
1005 // "Small": For |x| <= 1, approximate atan(x) directly by a polynomial
1006 // calculated using Sollya.
1007
1008 const Packet abs_x = pabs(x_in);
1009 const Packet x_signmask = pand(x_in, cst_signmask);
1010 const Packet large_mask = pcmp_lt(cst_one, abs_x);
1011 const Packet x = pselect(large_mask, preciprocal(abs_x), abs_x);
1012 const Packet p = patan_reduced_float(x);
1013 // Apply transformations according to the range reduction masks.
1014 Packet result = pselect(large_mask, psub(cst_pi_over_two, p), p);
1015 // Return correct sign
1016 return pxor(result, x_signmask);
1017}
1018
1019// Computes elementwise atan(x) for x in [-tan(pi/8):tan(pi/8)]
1020// with 2 ulp accuracy.
1021template <typename Packet>
1022EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet patan_reduced_double(const Packet& x) {
1023 const Packet q0 = pset1<Packet>(-0.33333333333330028569463365784031338989734649658203);
1024 const Packet q2 = pset1<Packet>(0.199999999990664090177006073645316064357757568359375);
1025 const Packet q4 = pset1<Packet>(-0.142857141937123677255527809393242932856082916259766);
1026 const Packet q6 = pset1<Packet>(0.111111065991039953404495577160560060292482376098633);
1027 const Packet q8 = pset1<Packet>(-9.0907812986129224452902519715280504897236824035645e-2);
1028 const Packet q10 = pset1<Packet>(7.6900542950704739442180368769186316058039665222168e-2);
1029 const Packet q12 = pset1<Packet>(-6.6410112986494976294871150912513257935643196105957e-2);
1030 const Packet q14 = pset1<Packet>(5.6920144995467943094258345126945641823112964630127e-2);
1031 const Packet q16 = pset1<Packet>(-4.3577020814990513608577771265117917209863662719727e-2);
1032 const Packet q18 = pset1<Packet>(2.1244050233624342527427586446719942614436149597168e-2);
1033
1034 // Approximate atan(x) on [0:tan(pi/8)] by a polynomial of the form
1035 // P(x) = x + x^3 * Q(x^2),
1036 // where Q(x^2) is a 9th order polynomial in x^2.
1037 // We evaluate even and odd terms in x^2 in parallel
1038 // to take advantage of instruction level parallelism
1039 // and hardware with multiple FMA units.
1040 const Packet x2 = pmul(x, x);
1041 const Packet x4 = pmul(x2, x2);
1042 Packet q_odd = pmadd(q18, x4, q14);
1043 Packet q_even = pmadd(q16, x4, q12);
1044 q_odd = pmadd(q_odd, x4, q10);
1045 q_even = pmadd(q_even, x4, q8);
1046 q_odd = pmadd(q_odd, x4, q6);
1047 q_even = pmadd(q_even, x4, q4);
1048 q_odd = pmadd(q_odd, x4, q2);
1049 q_even = pmadd(q_even, x4, q0);
1050 const Packet p = pmadd(q_odd, x2, q_even);
1051 return pmadd(p, pmul(x, x2), x);
1052}
1053
1054template <typename Packet>
1055EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet patan_double(const Packet& x_in) {
1056 typedef typename unpacket_traits<Packet>::type Scalar;
1057 static_assert(std::is_same<Scalar, double>::value, "Scalar type must be double");
1058
1059 constexpr double kPiOverTwo = static_cast<double>(EIGEN_PI / 2);
1060 constexpr double kPiOverFour = static_cast<double>(EIGEN_PI / 4);
1061 constexpr double kTanPiOverEight = 0.4142135623730950488016887;
1062 constexpr double kTan3PiOverEight = 2.4142135623730950488016887;
1063
1064 const Packet cst_signmask = pset1<Packet>(-0.0);
1065 const Packet cst_one = pset1<Packet>(1.0);
1066 const Packet cst_pi_over_two = pset1<Packet>(kPiOverTwo);
1067 const Packet cst_pi_over_four = pset1<Packet>(kPiOverFour);
1068 const Packet cst_large = pset1<Packet>(kTan3PiOverEight);
1069 const Packet cst_medium = pset1<Packet>(kTanPiOverEight);
1070
1071 // Use the same range reduction strategy (to [0:tan(pi/8)]) as the
1072 // Cephes library:
1073 // "Large": For x >= tan(3*pi/8), use atan(1/x) = pi/2 - atan(x).
1074 // "Medium": For x in [tan(pi/8) : tan(3*pi/8)),
1075 // use atan(x) = pi/4 + atan((x-1)/(x+1)).
1076 // "Small": For x < tan(pi/8), approximate atan(x) directly by a polynomial
1077 // calculated using Sollya.
1078
1079 const Packet abs_x = pabs(x_in);
1080 const Packet x_signmask = pand(x_in, cst_signmask);
1081 const Packet large_mask = pcmp_lt(cst_large, abs_x);
1082 const Packet medium_mask = pandnot(pcmp_lt(cst_medium, abs_x), large_mask);
1083
1084 Packet x = abs_x;
1085 x = pselect(large_mask, preciprocal(abs_x), x);
1086 x = pselect(medium_mask, pdiv(psub(abs_x, cst_one), padd(abs_x, cst_one)), x);
1087
1088 // Compute approximation of p ~= atan(x') where x' is the argument reduced to
1089 // [0:tan(pi/8)].
1090 Packet p = patan_reduced_double(x);
1091
1092 // Apply transformations according to the range reduction masks.
1093 p = pselect(large_mask, psub(cst_pi_over_two, p), p);
1094 p = pselect(medium_mask, padd(cst_pi_over_four, p), p);
1095 // Return the correct sign
1096 return pxor(p, x_signmask);
1097}
1098
1109template <typename T>
1110EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS T ptanh_float(const T& a_x) {
1111 // Clamp the inputs to the range [-c, c]
1112#ifdef EIGEN_VECTORIZE_FMA
1113 const T plus_clamp = pset1<T>(7.99881172180175781f);
1114 const T minus_clamp = pset1<T>(-7.99881172180175781f);
1115#else
1116 const T plus_clamp = pset1<T>(7.90531110763549805f);
1117 const T minus_clamp = pset1<T>(-7.90531110763549805f);
1118#endif
1119 const T tiny = pset1<T>(0.0004f);
1120 const T x = pmax(pmin(a_x, plus_clamp), minus_clamp);
1121 const T tiny_mask = pcmp_lt(pabs(a_x), tiny);
1122 // The monomial coefficients of the numerator polynomial (odd).
1123 const T alpha_1 = pset1<T>(4.89352455891786e-03f);
1124 const T alpha_3 = pset1<T>(6.37261928875436e-04f);
1125 const T alpha_5 = pset1<T>(1.48572235717979e-05f);
1126 const T alpha_7 = pset1<T>(5.12229709037114e-08f);
1127 const T alpha_9 = pset1<T>(-8.60467152213735e-11f);
1128 const T alpha_11 = pset1<T>(2.00018790482477e-13f);
1129 const T alpha_13 = pset1<T>(-2.76076847742355e-16f);
1130
1131 // The monomial coefficients of the denominator polynomial (even).
1132 const T beta_0 = pset1<T>(4.89352518554385e-03f);
1133 const T beta_2 = pset1<T>(2.26843463243900e-03f);
1134 const T beta_4 = pset1<T>(1.18534705686654e-04f);
1135 const T beta_6 = pset1<T>(1.19825839466702e-06f);
1136
1137 // Since the polynomials are odd/even, we need x^2.
1138 const T x2 = pmul(x, x);
1139
1140 // Evaluate the numerator polynomial p.
1141 T p = pmadd(x2, alpha_13, alpha_11);
1142 p = pmadd(x2, p, alpha_9);
1143 p = pmadd(x2, p, alpha_7);
1144 p = pmadd(x2, p, alpha_5);
1145 p = pmadd(x2, p, alpha_3);
1146 p = pmadd(x2, p, alpha_1);
1147 p = pmul(x, p);
1148
1149 // Evaluate the denominator polynomial q.
1150 T q = pmadd(x2, beta_6, beta_4);
1151 q = pmadd(x2, q, beta_2);
1152 q = pmadd(x2, q, beta_0);
1153
1154 // Divide the numerator by the denominator.
1155 return pselect(tiny_mask, x, pdiv(p, q));
1156}
1157
1158template <typename Packet>
1159EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet patanh_float(const Packet& x) {
1160 typedef typename unpacket_traits<Packet>::type Scalar;
1161 static_assert(std::is_same<Scalar, float>::value, "Scalar type must be float");
1162 const Packet half = pset1<Packet>(0.5f);
1163 const Packet x_gt_half = pcmp_le(half, pabs(x));
1164 // For |x| in [0:0.5] we use a polynomial approximation of the form
1165 // P(x) = x + x^3*(c3 + x^2 * (c5 + x^2 * (... x^2 * c11) ... )).
1166 const Packet C3 = pset1<Packet>(0.3333373963832855224609375f);
1167 const Packet C5 = pset1<Packet>(0.1997792422771453857421875f);
1168 const Packet C7 = pset1<Packet>(0.14672131836414337158203125f);
1169 const Packet C9 = pset1<Packet>(8.2311116158962249755859375e-2f);
1170 const Packet C11 = pset1<Packet>(0.1819281280040740966796875f);
1171 const Packet x2 = pmul(x, x);
1172 Packet p = pmadd(C11, x2, C9);
1173 p = pmadd(x2, p, C7);
1174 p = pmadd(x2, p, C5);
1175 p = pmadd(x2, p, C3);
1176 p = pmadd(pmul(x, x2), p, x);
1177
1178 // For |x| in ]0.5:1.0] we use atanh = 0.5*ln((1+x)/(1-x));
1179 const Packet one = pset1<Packet>(1.0f);
1180 Packet r = pdiv(padd(one, x), psub(one, x));
1181 r = pmul(half, plog(r));
1182 return pselect(x_gt_half, r, p);
1183}
1184
1185template <typename Packet>
1186EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet pdiv_complex(const Packet& x, const Packet& y) {
1187 typedef typename unpacket_traits<Packet>::as_real RealPacket;
1188 // In the following we annotate the code for the case where the inputs
1189 // are a pair length-2 SIMD vectors representing a single pair of complex
1190 // numbers x = a + i*b, y = c + i*d.
1191 const RealPacket y_abs = pabs(y.v); // |c|, |d|
1192 const RealPacket y_abs_flip = pcplxflip(Packet(y_abs)).v; // |d|, |c|
1193 const RealPacket y_max = pmax(y_abs, y_abs_flip); // max(|c|, |d|), max(|c|, |d|)
1194 const RealPacket y_scaled = pdiv(y.v, y_max); // c / max(|c|, |d|), d / max(|c|, |d|)
1195 // Compute scaled denominator.
1196 const RealPacket y_scaled_sq = pmul(y_scaled, y_scaled); // c'**2, d'**2
1197 const RealPacket denom = padd(y_scaled_sq, pcplxflip(Packet(y_scaled_sq)).v);
1198 Packet result_scaled = pmul(x, pconj(Packet(y_scaled))); // a * c' + b * d', -a * d + b * c
1199 // Divide elementwise by denom.
1200 result_scaled = Packet(pdiv(result_scaled.v, denom));
1201 // Rescale result
1202 return Packet(pdiv(result_scaled.v, y_max));
1203}
1204
1205template <typename Packet>
1206EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet plog_complex(const Packet& x) {
1207 typedef typename unpacket_traits<Packet>::type Scalar;
1208 typedef typename Scalar::value_type RealScalar;
1209 typedef typename unpacket_traits<Packet>::as_real RealPacket;
1210
1211 RealPacket real_mask_rp = peven_mask(x.v);
1212 Packet real_mask(real_mask_rp);
1213
1214 // Real part
1215 RealPacket x_flip = pcplxflip(x).v; // b, a
1216 Packet x_norm = phypot_complex(x); // sqrt(a^2 + b^2), sqrt(a^2 + b^2)
1217 RealPacket xlogr = plog(x_norm.v); // log(sqrt(a^2 + b^2)), log(sqrt(a^2 + b^2))
1218
1219 // Imag part
1220 RealPacket ximg = patan2(x.v, x_flip); // atan2(a, b), atan2(b, a)
1221
1222 const RealPacket cst_pos_inf = pset1<RealPacket>(NumTraits<RealScalar>::infinity());
1223 RealPacket x_abs = pabs(x.v);
1224 RealPacket is_x_pos_inf = pcmp_eq(x_abs, cst_pos_inf);
1225 RealPacket is_y_pos_inf = pcplxflip(Packet(is_x_pos_inf)).v;
1226 RealPacket is_any_inf = por(is_x_pos_inf, is_y_pos_inf);
1227 RealPacket xreal = pselect(is_any_inf, cst_pos_inf, xlogr);
1228
1229 Packet xres = pselect(real_mask, Packet(xreal), Packet(ximg)); // log(sqrt(a^2 + b^2)), atan2(b, a)
1230 return xres;
1231}
1232
1233template <typename Packet>
1234EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet pexp_complex(const Packet& a) {
1235 typedef typename unpacket_traits<Packet>::as_real RealPacket;
1236 typedef typename unpacket_traits<Packet>::type Scalar;
1237 typedef typename Scalar::value_type RealScalar;
1238 const RealPacket even_mask = peven_mask(a.v);
1239 const RealPacket odd_mask = pcplxflip(Packet(even_mask)).v;
1240
1241 // Let a = x + iy.
1242 // exp(a) = exp(x) * cis(y), plus some special edge-case handling.
1243
1244 // exp(x):
1245 RealPacket x = pand(a.v, even_mask);
1246 x = por(x, pcplxflip(Packet(x)).v);
1247 RealPacket expx = pexp(x); // exp(x);
1248
1249 // cis(y):
1250 RealPacket y = pand(odd_mask, a.v);
1251 y = por(y, pcplxflip(Packet(y)).v);
1252 RealPacket cisy = psincos_float<false, RealPacket, true>(y);
1253 cisy = pcplxflip(Packet(cisy)).v; // cos(y) + i * sin(y)
1254
1255 const RealPacket cst_pos_inf = pset1<RealPacket>(NumTraits<RealScalar>::infinity());
1256 const RealPacket cst_neg_inf = pset1<RealPacket>(-NumTraits<RealScalar>::infinity());
1257
1258 // If x is -inf, we know that cossin(y) is bounded,
1259 // so the result is (0, +/-0), where the sign of the imaginary part comes
1260 // from the sign of cossin(y).
1261 RealPacket cisy_sign = por(pandnot(cisy, pabs(cisy)), pset1<RealPacket>(RealScalar(1)));
1262 cisy = pselect(pcmp_eq(x, cst_neg_inf), cisy_sign, cisy);
1263
1264 // If x is inf, and cos(y) has unknown sign (y is inf or NaN), the result
1265 // is (+/-inf, NaN), where the signs are undetermined (take the sign of y).
1266 RealPacket y_sign = por(pandnot(y, pabs(y)), pset1<RealPacket>(RealScalar(1)));
1267 cisy = pselect(pand(pcmp_eq(x, cst_pos_inf), pisnan(cisy)), pand(y_sign, even_mask), cisy);
1268 Packet result = Packet(pmul(expx, cisy));
1269
1270 // If y is +/- 0, the input is real, so take the real result for consistency.
1271 result = pselect(Packet(pcmp_eq(y, pzero(y))), Packet(por(pand(expx, even_mask), pand(y, odd_mask))), result);
1272
1273 return result;
1274}
1275
1276template <typename Packet>
1277EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet psqrt_complex(const Packet& a) {
1278 typedef typename unpacket_traits<Packet>::type Scalar;
1279 typedef typename Scalar::value_type RealScalar;
1280 typedef typename unpacket_traits<Packet>::as_real RealPacket;
1281
1282 // Computes the principal sqrt of the complex numbers in the input.
1283 //
1284 // For example, for packets containing 2 complex numbers stored in interleaved format
1285 // a = [a0, a1] = [x0, y0, x1, y1],
1286 // where x0 = real(a0), y0 = imag(a0) etc., this function returns
1287 // b = [b0, b1] = [u0, v0, u1, v1],
1288 // such that b0^2 = a0, b1^2 = a1.
1289 //
1290 // To derive the formula for the complex square roots, let's consider the equation for
1291 // a single complex square root of the number x + i*y. We want to find real numbers
1292 // u and v such that
1293 // (u + i*v)^2 = x + i*y <=>
1294 // u^2 - v^2 + i*2*u*v = x + i*v.
1295 // By equating the real and imaginary parts we get:
1296 // u^2 - v^2 = x
1297 // 2*u*v = y.
1298 //
1299 // For x >= 0, this has the numerically stable solution
1300 // u = sqrt(0.5 * (x + sqrt(x^2 + y^2)))
1301 // v = 0.5 * (y / u)
1302 // and for x < 0,
1303 // v = sign(y) * sqrt(0.5 * (-x + sqrt(x^2 + y^2)))
1304 // u = 0.5 * (y / v)
1305 //
1306 // To avoid unnecessary over- and underflow, we compute sqrt(x^2 + y^2) as
1307 // l = max(|x|, |y|) * sqrt(1 + (min(|x|, |y|) / max(|x|, |y|))^2) ,
1308
1309 // In the following, without lack of generality, we have annotated the code, assuming
1310 // that the input is a packet of 2 complex numbers.
1311 //
1312 // Step 1. Compute l = [l0, l0, l1, l1], where
1313 // l0 = sqrt(x0^2 + y0^2), l1 = sqrt(x1^2 + y1^2)
1314 // To avoid over- and underflow, we use the stable formula for each hypotenuse
1315 // l0 = (min0 == 0 ? max0 : max0 * sqrt(1 + (min0/max0)**2)),
1316 // where max0 = max(|x0|, |y0|), min0 = min(|x0|, |y0|), and similarly for l1.
1317
1318 RealPacket a_abs = pabs(a.v); // [|x0|, |y0|, |x1|, |y1|]
1319 RealPacket a_abs_flip = pcplxflip(Packet(a_abs)).v; // [|y0|, |x0|, |y1|, |x1|]
1320 RealPacket a_max = pmax(a_abs, a_abs_flip);
1321 RealPacket a_min = pmin(a_abs, a_abs_flip);
1322 RealPacket a_min_zero_mask = pcmp_eq(a_min, pzero(a_min));
1323 RealPacket a_max_zero_mask = pcmp_eq(a_max, pzero(a_max));
1324 RealPacket r = pdiv(a_min, a_max);
1325 const RealPacket cst_one = pset1<RealPacket>(RealScalar(1));
1326 RealPacket l = pmul(a_max, psqrt(padd(cst_one, pmul(r, r)))); // [l0, l0, l1, l1]
1327 // Set l to a_max if a_min is zero.
1328 l = pselect(a_min_zero_mask, a_max, l);
1329
1330 // Step 2. Compute [rho0, *, rho1, *], where
1331 // rho0 = sqrt(0.5 * (l0 + |x0|)), rho1 = sqrt(0.5 * (l1 + |x1|))
1332 // We don't care about the imaginary parts computed here. They will be overwritten later.
1333 const RealPacket cst_half = pset1<RealPacket>(RealScalar(0.5));
1334 Packet rho;
1335 rho.v = psqrt(pmul(cst_half, padd(a_abs, l)));
1336
1337 // Step 3. Compute [rho0, eta0, rho1, eta1], where
1338 // eta0 = (y0 / l0) / 2, and eta1 = (y1 / l1) / 2.
1339 // set eta = 0 of input is 0 + i0.
1340 RealPacket eta = pandnot(pmul(cst_half, pdiv(a.v, pcplxflip(rho).v)), a_max_zero_mask);
1341 RealPacket real_mask = peven_mask(a.v);
1342 Packet positive_real_result;
1343 // Compute result for inputs with positive real part.
1344 positive_real_result.v = pselect(real_mask, rho.v, eta);
1345
1346 // Step 4. Compute solution for inputs with negative real part:
1347 // [|eta0|, sign(y0)*rho0, |eta1|, sign(y1)*rho1]
1348 const RealPacket cst_imag_sign_mask = pset1<Packet>(Scalar(RealScalar(0.0), RealScalar(-0.0))).v;
1349 RealPacket imag_signs = pand(a.v, cst_imag_sign_mask);
1350 Packet negative_real_result;
1351 // Notice that rho is positive, so taking it's absolute value is a noop.
1352 negative_real_result.v = por(pabs(pcplxflip(positive_real_result).v), imag_signs);
1353
1354 // Step 5. Select solution branch based on the sign of the real parts.
1355 Packet negative_real_mask;
1356 negative_real_mask.v = pcmp_lt(pand(real_mask, a.v), pzero(a.v));
1357 negative_real_mask.v = por(negative_real_mask.v, pcplxflip(negative_real_mask).v);
1358 Packet result = pselect(negative_real_mask, negative_real_result, positive_real_result);
1359
1360 // Step 6. Handle special cases for infinities:
1361 // * If z is (x,+∞), the result is (+∞,+∞) even if x is NaN
1362 // * If z is (x,-∞), the result is (+∞,-∞) even if x is NaN
1363 // * If z is (-∞,y), the result is (0*|y|,+∞) for finite or NaN y
1364 // * If z is (+∞,y), the result is (+∞,0*|y|) for finite or NaN y
1365 const RealPacket cst_pos_inf = pset1<RealPacket>(NumTraits<RealScalar>::infinity());
1366 Packet is_inf;
1367 is_inf.v = pcmp_eq(a_abs, cst_pos_inf);
1368 Packet is_real_inf;
1369 is_real_inf.v = pand(is_inf.v, real_mask);
1370 is_real_inf = por(is_real_inf, pcplxflip(is_real_inf));
1371 // prepare packet of (+∞,0*|y|) or (0*|y|,+∞), depending on the sign of the infinite real part.
1372 Packet real_inf_result;
1373 real_inf_result.v = pmul(a_abs, pset1<Packet>(Scalar(RealScalar(1.0), RealScalar(0.0))).v);
1374 real_inf_result.v = pselect(negative_real_mask.v, pcplxflip(real_inf_result).v, real_inf_result.v);
1375 // prepare packet of (+∞,+∞) or (+∞,-∞), depending on the sign of the infinite imaginary part.
1376 Packet is_imag_inf;
1377 is_imag_inf.v = pandnot(is_inf.v, real_mask);
1378 is_imag_inf = por(is_imag_inf, pcplxflip(is_imag_inf));
1379 Packet imag_inf_result;
1380 imag_inf_result.v = por(pand(cst_pos_inf, real_mask), pandnot(a.v, real_mask));
1381 // unless otherwise specified, if either the real or imaginary component is nan, the entire result is nan
1382 Packet result_is_nan = pisnan(result);
1383 result = por(result_is_nan, result);
1384
1385 return pselect(is_imag_inf, imag_inf_result, pselect(is_real_inf, real_inf_result, result));
1386}
1387
1388// \internal \returns the norm of a complex number z = x + i*y, defined as sqrt(x^2 + y^2).
1389// Implemented using the hypot(a,b) algorithm from https://doi.org/10.48550/arXiv.1904.09481
1390template <typename Packet>
1391EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet phypot_complex(const Packet& a) {
1392 typedef typename unpacket_traits<Packet>::type Scalar;
1393 typedef typename Scalar::value_type RealScalar;
1394 typedef typename unpacket_traits<Packet>::as_real RealPacket;
1395
1396 const RealPacket cst_zero_rp = pset1<RealPacket>(static_cast<RealScalar>(0.0));
1397 const RealPacket cst_minus_one_rp = pset1<RealPacket>(static_cast<RealScalar>(-1.0));
1398 const RealPacket cst_two_rp = pset1<RealPacket>(static_cast<RealScalar>(2.0));
1399 const RealPacket evenmask = peven_mask(a.v);
1400
1401 RealPacket a_abs = pabs(a.v);
1402 RealPacket a_flip = pcplxflip(Packet(a_abs)).v; // |b|, |a|
1403 RealPacket a_all = pselect(evenmask, a_abs, a_flip); // |a|, |a|
1404 RealPacket b_all = pselect(evenmask, a_flip, a_abs); // |b|, |b|
1405
1406 RealPacket a2 = pmul(a.v, a.v); // |a^2, b^2|
1407 RealPacket a2_flip = pcplxflip(Packet(a2)).v; // |b^2, a^2|
1408 RealPacket h = psqrt(padd(a2, a2_flip)); // |sqrt(a^2 + b^2), sqrt(a^2 + b^2)|
1409 RealPacket h_sq = pmul(h, h); // |a^2 + b^2, a^2 + b^2|
1410 RealPacket a_sq = pselect(evenmask, a2, a2_flip); // |a^2, a^2|
1411 RealPacket m_h_sq = pmul(h_sq, cst_minus_one_rp);
1412 RealPacket m_a_sq = pmul(a_sq, cst_minus_one_rp);
1413 RealPacket x = psub(psub(pmadd(h, h, m_h_sq), pmadd(b_all, b_all, psub(a_sq, h_sq))), pmadd(a_all, a_all, m_a_sq));
1414 h = psub(h, pdiv(x, pmul(cst_two_rp, h))); // |h - x/(2*h), h - x/(2*h)|
1415
1416 // handle zero-case
1417 RealPacket iszero = pcmp_eq(por(a_abs, a_flip), cst_zero_rp);
1418
1419 h = pandnot(h, iszero); // |sqrt(a^2+b^2), sqrt(a^2+b^2)|
1420 return Packet(h); // |sqrt(a^2+b^2), sqrt(a^2+b^2)|
1421}
1422
1423template <typename Packet>
1424struct psign_impl<Packet, std::enable_if_t<!NumTraits<typename unpacket_traits<Packet>::type>::IsComplex &&
1425 !NumTraits<typename unpacket_traits<Packet>::type>::IsInteger>> {
1426 static EIGEN_DEVICE_FUNC inline Packet run(const Packet& a) {
1427 using Scalar = typename unpacket_traits<Packet>::type;
1428 const Packet cst_one = pset1<Packet>(Scalar(1));
1429 const Packet cst_zero = pzero(a);
1430
1431 const Packet abs_a = pabs(a);
1432 const Packet sign_mask = pandnot(a, abs_a);
1433 const Packet nonzero_mask = pcmp_lt(cst_zero, abs_a);
1434
1435 return pselect(nonzero_mask, por(sign_mask, cst_one), abs_a);
1436 }
1437};
1438
1439template <typename Packet>
1440struct psign_impl<Packet, std::enable_if_t<!NumTraits<typename unpacket_traits<Packet>::type>::IsComplex &&
1441 NumTraits<typename unpacket_traits<Packet>::type>::IsSigned &&
1442 NumTraits<typename unpacket_traits<Packet>::type>::IsInteger>> {
1443 static EIGEN_DEVICE_FUNC inline Packet run(const Packet& a) {
1444 using Scalar = typename unpacket_traits<Packet>::type;
1445 const Packet cst_one = pset1<Packet>(Scalar(1));
1446 const Packet cst_minus_one = pset1<Packet>(Scalar(-1));
1447 const Packet cst_zero = pzero(a);
1448
1449 const Packet positive_mask = pcmp_lt(cst_zero, a);
1450 const Packet positive = pand(positive_mask, cst_one);
1451 const Packet negative_mask = pcmp_lt(a, cst_zero);
1452 const Packet negative = pand(negative_mask, cst_minus_one);
1453
1454 return por(positive, negative);
1455 }
1456};
1457
1458template <typename Packet>
1459struct psign_impl<Packet, std::enable_if_t<!NumTraits<typename unpacket_traits<Packet>::type>::IsComplex &&
1460 !NumTraits<typename unpacket_traits<Packet>::type>::IsSigned &&
1461 NumTraits<typename unpacket_traits<Packet>::type>::IsInteger>> {
1462 static EIGEN_DEVICE_FUNC inline Packet run(const Packet& a) {
1463 using Scalar = typename unpacket_traits<Packet>::type;
1464 const Packet cst_one = pset1<Packet>(Scalar(1));
1465 const Packet cst_zero = pzero(a);
1466
1467 const Packet zero_mask = pcmp_eq(cst_zero, a);
1468 return pandnot(cst_one, zero_mask);
1469 }
1470};
1471
1472// \internal \returns the the sign of a complex number z, defined as z / abs(z).
1473template <typename Packet>
1474struct psign_impl<Packet, std::enable_if_t<NumTraits<typename unpacket_traits<Packet>::type>::IsComplex &&
1475 unpacket_traits<Packet>::vectorizable>> {
1476 static EIGEN_DEVICE_FUNC inline Packet run(const Packet& a) {
1477 typedef typename unpacket_traits<Packet>::type Scalar;
1478 typedef typename Scalar::value_type RealScalar;
1479 typedef typename unpacket_traits<Packet>::as_real RealPacket;
1480
1481 // Step 1. Compute (for each element z = x + i*y in a)
1482 // l = abs(z) = sqrt(x^2 + y^2).
1483 // To avoid over- and underflow, we use the stable formula for each hypotenuse
1484 // l = (zmin == 0 ? zmax : zmax * sqrt(1 + (zmin/zmax)**2)),
1485 // where zmax = max(|x|, |y|), zmin = min(|x|, |y|),
1486 RealPacket a_abs = pabs(a.v);
1487 RealPacket a_abs_flip = pcplxflip(Packet(a_abs)).v;
1488 RealPacket a_max = pmax(a_abs, a_abs_flip);
1489 RealPacket a_min = pmin(a_abs, a_abs_flip);
1490 RealPacket a_min_zero_mask = pcmp_eq(a_min, pzero(a_min));
1491 RealPacket a_max_zero_mask = pcmp_eq(a_max, pzero(a_max));
1492 RealPacket r = pdiv(a_min, a_max);
1493 const RealPacket cst_one = pset1<RealPacket>(RealScalar(1));
1494 RealPacket l = pmul(a_max, psqrt(padd(cst_one, pmul(r, r)))); // [l0, l0, l1, l1]
1495 // Set l to a_max if a_min is zero, since the roundtrip sqrt(a_max^2) may be
1496 // lossy.
1497 l = pselect(a_min_zero_mask, a_max, l);
1498 // Step 2 compute a / abs(a).
1499 RealPacket sign_as_real = pandnot(pdiv(a.v, l), a_max_zero_mask);
1500 Packet sign;
1501 sign.v = sign_as_real;
1502 return sign;
1503 }
1504};
1505
1506// TODO(rmlarsen): The following set of utilities for double word arithmetic
1507// should perhaps be refactored as a separate file, since it would be generally
1508// useful for special function implementation etc. Writing the algorithms in
1509// terms if a double word type would also make the code more readable.
1510
1511// This function splits x into the nearest integer n and fractional part r,
1512// such that x = n + r holds exactly.
1513template <typename Packet>
1514EIGEN_STRONG_INLINE void absolute_split(const Packet& x, Packet& n, Packet& r) {
1515 n = pround(x);
1516 r = psub(x, n);
1517}
1518
1519// This function computes the sum {s, r}, such that x + y = s_hi + s_lo
1520// holds exactly, and s_hi = fl(x+y), if |x| >= |y|.
1521template <typename Packet>
1522EIGEN_STRONG_INLINE void fast_twosum(const Packet& x, const Packet& y, Packet& s_hi, Packet& s_lo) {
1523 s_hi = padd(x, y);
1524 const Packet t = psub(s_hi, x);
1525 s_lo = psub(y, t);
1526}
1527
1528#ifdef EIGEN_VECTORIZE_FMA
1529// This function implements the extended precision product of
1530// a pair of floating point numbers. Given {x, y}, it computes the pair
1531// {p_hi, p_lo} such that x * y = p_hi + p_lo holds exactly and
1532// p_hi = fl(x * y).
1533template <typename Packet>
1534EIGEN_STRONG_INLINE void twoprod(const Packet& x, const Packet& y, Packet& p_hi, Packet& p_lo) {
1535 p_hi = pmul(x, y);
1536 p_lo = pmsub(x, y, p_hi);
1537}
1538
1539#else
1540
1541// This function implements the Veltkamp splitting. Given a floating point
1542// number x it returns the pair {x_hi, x_lo} such that x_hi + x_lo = x holds
1543// exactly and that half of the significant of x fits in x_hi.
1544// This is Algorithm 3 from Jean-Michel Muller, "Elementary Functions",
1545// 3rd edition, Birkh\"auser, 2016.
1546template <typename Packet>
1547EIGEN_STRONG_INLINE void veltkamp_splitting(const Packet& x, Packet& x_hi, Packet& x_lo) {
1548 typedef typename unpacket_traits<Packet>::type Scalar;
1549 EIGEN_CONSTEXPR int shift = (NumTraits<Scalar>::digits() + 1) / 2;
1550 const Scalar shift_scale = Scalar(uint64_t(1) << shift); // Scalar constructor not necessarily constexpr.
1551 const Packet gamma = pmul(pset1<Packet>(shift_scale + Scalar(1)), x);
1552 Packet rho = psub(x, gamma);
1553 x_hi = padd(rho, gamma);
1554 x_lo = psub(x, x_hi);
1555}
1556
1557// This function implements Dekker's algorithm for products x * y.
1558// Given floating point numbers {x, y} computes the pair
1559// {p_hi, p_lo} such that x * y = p_hi + p_lo holds exactly and
1560// p_hi = fl(x * y).
1561template <typename Packet>
1562EIGEN_STRONG_INLINE void twoprod(const Packet& x, const Packet& y, Packet& p_hi, Packet& p_lo) {
1563 Packet x_hi, x_lo, y_hi, y_lo;
1564 veltkamp_splitting(x, x_hi, x_lo);
1565 veltkamp_splitting(y, y_hi, y_lo);
1566
1567 p_hi = pmul(x, y);
1568 p_lo = pmadd(x_hi, y_hi, pnegate(p_hi));
1569 p_lo = pmadd(x_hi, y_lo, p_lo);
1570 p_lo = pmadd(x_lo, y_hi, p_lo);
1571 p_lo = pmadd(x_lo, y_lo, p_lo);
1572}
1573
1574#endif // EIGEN_VECTORIZE_FMA
1575
1576// This function implements Dekker's algorithm for the addition
1577// of two double word numbers represented by {x_hi, x_lo} and {y_hi, y_lo}.
1578// It returns the result as a pair {s_hi, s_lo} such that
1579// x_hi + x_lo + y_hi + y_lo = s_hi + s_lo holds exactly.
1580// This is Algorithm 5 from Jean-Michel Muller, "Elementary Functions",
1581// 3rd edition, Birkh\"auser, 2016.
1582template <typename Packet>
1583EIGEN_STRONG_INLINE void twosum(const Packet& x_hi, const Packet& x_lo, const Packet& y_hi, const Packet& y_lo,
1584 Packet& s_hi, Packet& s_lo) {
1585 const Packet x_greater_mask = pcmp_lt(pabs(y_hi), pabs(x_hi));
1586 Packet r_hi_1, r_lo_1;
1587 fast_twosum(x_hi, y_hi, r_hi_1, r_lo_1);
1588 Packet r_hi_2, r_lo_2;
1589 fast_twosum(y_hi, x_hi, r_hi_2, r_lo_2);
1590 const Packet r_hi = pselect(x_greater_mask, r_hi_1, r_hi_2);
1591
1592 const Packet s1 = padd(padd(y_lo, r_lo_1), x_lo);
1593 const Packet s2 = padd(padd(x_lo, r_lo_2), y_lo);
1594 const Packet s = pselect(x_greater_mask, s1, s2);
1595
1596 fast_twosum(r_hi, s, s_hi, s_lo);
1597}
1598
1599// This is a version of twosum for double word numbers,
1600// which assumes that |x_hi| >= |y_hi|.
1601template <typename Packet>
1602EIGEN_STRONG_INLINE void fast_twosum(const Packet& x_hi, const Packet& x_lo, const Packet& y_hi, const Packet& y_lo,
1603 Packet& s_hi, Packet& s_lo) {
1604 Packet r_hi, r_lo;
1605 fast_twosum(x_hi, y_hi, r_hi, r_lo);
1606 const Packet s = padd(padd(y_lo, r_lo), x_lo);
1607 fast_twosum(r_hi, s, s_hi, s_lo);
1608}
1609
1610// This is a version of twosum for adding a floating point number x to
1611// double word number {y_hi, y_lo} number, with the assumption
1612// that |x| >= |y_hi|.
1613template <typename Packet>
1614EIGEN_STRONG_INLINE void fast_twosum(const Packet& x, const Packet& y_hi, const Packet& y_lo, Packet& s_hi,
1615 Packet& s_lo) {
1616 Packet r_hi, r_lo;
1617 fast_twosum(x, y_hi, r_hi, r_lo);
1618 const Packet s = padd(y_lo, r_lo);
1619 fast_twosum(r_hi, s, s_hi, s_lo);
1620}
1621
1622// This function implements the multiplication of a double word
1623// number represented by {x_hi, x_lo} by a floating point number y.
1624// It returns the result as a pair {p_hi, p_lo} such that
1625// (x_hi + x_lo) * y = p_hi + p_lo hold with a relative error
1626// of less than 2*2^{-2p}, where p is the number of significand bit
1627// in the floating point type.
1628// This is Algorithm 7 from Jean-Michel Muller, "Elementary Functions",
1629// 3rd edition, Birkh\"auser, 2016.
1630template <typename Packet>
1631EIGEN_STRONG_INLINE void twoprod(const Packet& x_hi, const Packet& x_lo, const Packet& y, Packet& p_hi, Packet& p_lo) {
1632 Packet c_hi, c_lo1;
1633 twoprod(x_hi, y, c_hi, c_lo1);
1634 const Packet c_lo2 = pmul(x_lo, y);
1635 Packet t_hi, t_lo1;
1636 fast_twosum(c_hi, c_lo2, t_hi, t_lo1);
1637 const Packet t_lo2 = padd(t_lo1, c_lo1);
1638 fast_twosum(t_hi, t_lo2, p_hi, p_lo);
1639}
1640
1641// This function implements the multiplication of two double word
1642// numbers represented by {x_hi, x_lo} and {y_hi, y_lo}.
1643// It returns the result as a pair {p_hi, p_lo} such that
1644// (x_hi + x_lo) * (y_hi + y_lo) = p_hi + p_lo holds with a relative error
1645// of less than 2*2^{-2p}, where p is the number of significand bit
1646// in the floating point type.
1647template <typename Packet>
1648EIGEN_STRONG_INLINE void twoprod(const Packet& x_hi, const Packet& x_lo, const Packet& y_hi, const Packet& y_lo,
1649 Packet& p_hi, Packet& p_lo) {
1650 Packet p_hi_hi, p_hi_lo;
1651 twoprod(x_hi, x_lo, y_hi, p_hi_hi, p_hi_lo);
1652 Packet p_lo_hi, p_lo_lo;
1653 twoprod(x_hi, x_lo, y_lo, p_lo_hi, p_lo_lo);
1654 fast_twosum(p_hi_hi, p_hi_lo, p_lo_hi, p_lo_lo, p_hi, p_lo);
1655}
1656
1657// This function implements the division of double word {x_hi, x_lo}
1658// by float y. This is Algorithm 15 from "Tight and rigourous error bounds
1659// for basic building blocks of double-word arithmetic", Joldes, Muller, & Popescu,
1660// 2017. https://hal.archives-ouvertes.fr/hal-01351529
1661template <typename Packet>
1662void doubleword_div_fp(const Packet& x_hi, const Packet& x_lo, const Packet& y, Packet& z_hi, Packet& z_lo) {
1663 const Packet t_hi = pdiv(x_hi, y);
1664 Packet pi_hi, pi_lo;
1665 twoprod(t_hi, y, pi_hi, pi_lo);
1666 const Packet delta_hi = psub(x_hi, pi_hi);
1667 const Packet delta_t = psub(delta_hi, pi_lo);
1668 const Packet delta = padd(delta_t, x_lo);
1669 const Packet t_lo = pdiv(delta, y);
1670 fast_twosum(t_hi, t_lo, z_hi, z_lo);
1671}
1672
1673// This function computes log2(x) and returns the result as a double word.
1674template <typename Scalar>
1675struct accurate_log2 {
1676 template <typename Packet>
1677 EIGEN_STRONG_INLINE void operator()(const Packet& x, Packet& log2_x_hi, Packet& log2_x_lo) {
1678 log2_x_hi = plog2(x);
1679 log2_x_lo = pzero(x);
1680 }
1681};
1682
1683// This specialization uses a more accurate algorithm to compute log2(x) for
1684// floats in [1/sqrt(2);sqrt(2)] with a relative accuracy of ~6.42e-10.
1685// This additional accuracy is needed to counter the error-magnification
1686// inherent in multiplying by a potentially large exponent in pow(x,y).
1687// The minimax polynomial used was calculated using the Sollya tool.
1688// See sollya.org.
1689template <>
1690struct accurate_log2<float> {
1691 template <typename Packet>
1692 EIGEN_STRONG_INLINE void operator()(const Packet& z, Packet& log2_x_hi, Packet& log2_x_lo) {
1693 // The function log(1+x)/x is approximated in the interval
1694 // [1/sqrt(2)-1;sqrt(2)-1] by a degree 10 polynomial of the form
1695 // Q(x) = (C0 + x * (C1 + x * (C2 + x * (C3 + x * P(x))))),
1696 // where the degree 6 polynomial P(x) is evaluated in single precision,
1697 // while the remaining 4 terms of Q(x), as well as the final multiplication by x
1698 // to reconstruct log(1+x) are evaluated in extra precision using
1699 // double word arithmetic. C0 through C3 are extra precise constants
1700 // stored as double words.
1701 //
1702 // The polynomial coefficients were calculated using Sollya commands:
1703 // > n = 10;
1704 // > f = log2(1+x)/x;
1705 // > interval = [sqrt(0.5)-1;sqrt(2)-1];
1706 // > p = fpminimax(f,n,[|double,double,double,double,single...|],interval,relative,floating);
1707
1708 const Packet p6 = pset1<Packet>(9.703654795885e-2f);
1709 const Packet p5 = pset1<Packet>(-0.1690667718648f);
1710 const Packet p4 = pset1<Packet>(0.1720575392246f);
1711 const Packet p3 = pset1<Packet>(-0.1789081543684f);
1712 const Packet p2 = pset1<Packet>(0.2050433009862f);
1713 const Packet p1 = pset1<Packet>(-0.2404672354459f);
1714 const Packet p0 = pset1<Packet>(0.2885761857032f);
1715
1716 const Packet C3_hi = pset1<Packet>(-0.360674142838f);
1717 const Packet C3_lo = pset1<Packet>(-6.13283912543e-09f);
1718 const Packet C2_hi = pset1<Packet>(0.480897903442f);
1719 const Packet C2_lo = pset1<Packet>(-1.44861207474e-08f);
1720 const Packet C1_hi = pset1<Packet>(-0.721347510815f);
1721 const Packet C1_lo = pset1<Packet>(-4.84483164698e-09f);
1722 const Packet C0_hi = pset1<Packet>(1.44269502163f);
1723 const Packet C0_lo = pset1<Packet>(2.01711713999e-08f);
1724 const Packet one = pset1<Packet>(1.0f);
1725
1726 const Packet x = psub(z, one);
1727 // Evaluate P(x) in working precision.
1728 // We evaluate it in multiple parts to improve instruction level
1729 // parallelism.
1730 Packet x2 = pmul(x, x);
1731 Packet p_even = pmadd(p6, x2, p4);
1732 p_even = pmadd(p_even, x2, p2);
1733 p_even = pmadd(p_even, x2, p0);
1734 Packet p_odd = pmadd(p5, x2, p3);
1735 p_odd = pmadd(p_odd, x2, p1);
1736 Packet p = pmadd(p_odd, x, p_even);
1737
1738 // Now evaluate the low-order tems of Q(x) in double word precision.
1739 // In the following, due to the alternating signs and the fact that
1740 // |x| < sqrt(2)-1, we can assume that |C*_hi| >= q_i, and use
1741 // fast_twosum instead of the slower twosum.
1742 Packet q_hi, q_lo;
1743 Packet t_hi, t_lo;
1744 // C3 + x * p(x)
1745 twoprod(p, x, t_hi, t_lo);
1746 fast_twosum(C3_hi, C3_lo, t_hi, t_lo, q_hi, q_lo);
1747 // C2 + x * p(x)
1748 twoprod(q_hi, q_lo, x, t_hi, t_lo);
1749 fast_twosum(C2_hi, C2_lo, t_hi, t_lo, q_hi, q_lo);
1750 // C1 + x * p(x)
1751 twoprod(q_hi, q_lo, x, t_hi, t_lo);
1752 fast_twosum(C1_hi, C1_lo, t_hi, t_lo, q_hi, q_lo);
1753 // C0 + x * p(x)
1754 twoprod(q_hi, q_lo, x, t_hi, t_lo);
1755 fast_twosum(C0_hi, C0_lo, t_hi, t_lo, q_hi, q_lo);
1756
1757 // log(z) ~= x * Q(x)
1758 twoprod(q_hi, q_lo, x, log2_x_hi, log2_x_lo);
1759 }
1760};
1761
1762// This specialization uses a more accurate algorithm to compute log2(x) for
1763// floats in [1/sqrt(2);sqrt(2)] with a relative accuracy of ~1.27e-18.
1764// This additional accuracy is needed to counter the error-magnification
1765// inherent in multiplying by a potentially large exponent in pow(x,y).
1766// The minimax polynomial used was calculated using the Sollya tool.
1767// See sollya.org.
1768
1769template <>
1770struct accurate_log2<double> {
1771 template <typename Packet>
1772 EIGEN_STRONG_INLINE void operator()(const Packet& x, Packet& log2_x_hi, Packet& log2_x_lo) {
1773 // We use a transformation of variables:
1774 // r = c * (x-1) / (x+1),
1775 // such that
1776 // log2(x) = log2((1 + r/c) / (1 - r/c)) = f(r).
1777 // The function f(r) can be approximated well using an odd polynomial
1778 // of the form
1779 // P(r) = ((Q(r^2) * r^2 + C) * r^2 + 1) * r,
1780 // For the implementation of log2<double> here, Q is of degree 6 with
1781 // coefficient represented in working precision (double), while C is a
1782 // constant represented in extra precision as a double word to achieve
1783 // full accuracy.
1784 //
1785 // The polynomial coefficients were computed by the Sollya script:
1786 //
1787 // c = 2 / log(2);
1788 // trans = c * (x-1)/(x+1);
1789 // itrans = (1+x/c)/(1-x/c);
1790 // interval=[trans(sqrt(0.5)); trans(sqrt(2))];
1791 // print(interval);
1792 // f = log2(itrans(x));
1793 // p=fpminimax(f,[|1,3,5,7,9,11,13,15,17|],[|1,DD,double...|],interval,relative,floating);
1794 const Packet q12 = pset1<Packet>(2.87074255468000586e-9);
1795 const Packet q10 = pset1<Packet>(2.38957980901884082e-8);
1796 const Packet q8 = pset1<Packet>(2.31032094540014656e-7);
1797 const Packet q6 = pset1<Packet>(2.27279857398537278e-6);
1798 const Packet q4 = pset1<Packet>(2.31271023278625638e-5);
1799 const Packet q2 = pset1<Packet>(2.47556738444535513e-4);
1800 const Packet q0 = pset1<Packet>(2.88543873228900172e-3);
1801 const Packet C_hi = pset1<Packet>(0.0400377511598501157);
1802 const Packet C_lo = pset1<Packet>(-4.77726582251425391e-19);
1803 const Packet one = pset1<Packet>(1.0);
1804
1805 const Packet cst_2_log2e_hi = pset1<Packet>(2.88539008177792677);
1806 const Packet cst_2_log2e_lo = pset1<Packet>(4.07660016854549667e-17);
1807 // c * (x - 1)
1808 Packet t_hi, t_lo;
1809 // t = c * (x-1)
1810 twoprod(cst_2_log2e_hi, cst_2_log2e_lo, psub(x, one), t_hi, t_lo);
1811 // r = c * (x-1) / (x+1),
1812 Packet r_hi, r_lo;
1813 doubleword_div_fp(t_hi, t_lo, padd(x, one), r_hi, r_lo);
1814
1815 // r2 = r * r
1816 Packet r2_hi, r2_lo;
1817 twoprod(r_hi, r_lo, r_hi, r_lo, r2_hi, r2_lo);
1818 // r4 = r2 * r2
1819 Packet r4_hi, r4_lo;
1820 twoprod(r2_hi, r2_lo, r2_hi, r2_lo, r4_hi, r4_lo);
1821
1822 // Evaluate Q(r^2) in working precision. We evaluate it in two parts
1823 // (even and odd in r^2) to improve instruction level parallelism.
1824 Packet q_even = pmadd(q12, r4_hi, q8);
1825 Packet q_odd = pmadd(q10, r4_hi, q6);
1826 q_even = pmadd(q_even, r4_hi, q4);
1827 q_odd = pmadd(q_odd, r4_hi, q2);
1828 q_even = pmadd(q_even, r4_hi, q0);
1829 Packet q = pmadd(q_odd, r2_hi, q_even);
1830
1831 // Now evaluate the low order terms of P(x) in double word precision.
1832 // In the following, due to the increasing magnitude of the coefficients
1833 // and r being constrained to [-0.5, 0.5] we can use fast_twosum instead
1834 // of the slower twosum.
1835 // Q(r^2) * r^2
1836 Packet p_hi, p_lo;
1837 twoprod(r2_hi, r2_lo, q, p_hi, p_lo);
1838 // Q(r^2) * r^2 + C
1839 Packet p1_hi, p1_lo;
1840 fast_twosum(C_hi, C_lo, p_hi, p_lo, p1_hi, p1_lo);
1841 // (Q(r^2) * r^2 + C) * r^2
1842 Packet p2_hi, p2_lo;
1843 twoprod(r2_hi, r2_lo, p1_hi, p1_lo, p2_hi, p2_lo);
1844 // ((Q(r^2) * r^2 + C) * r^2 + 1)
1845 Packet p3_hi, p3_lo;
1846 fast_twosum(one, p2_hi, p2_lo, p3_hi, p3_lo);
1847
1848 // log(z) ~= ((Q(r^2) * r^2 + C) * r^2 + 1) * r
1849 twoprod(p3_hi, p3_lo, r_hi, r_lo, log2_x_hi, log2_x_lo);
1850 }
1851};
1852
1853// This function computes exp2(x) (i.e. 2**x).
1854template <typename Scalar>
1855struct fast_accurate_exp2 {
1856 template <typename Packet>
1857 EIGEN_STRONG_INLINE Packet operator()(const Packet& x) {
1858 // TODO(rmlarsen): Add a pexp2 packetop.
1859 return pexp(pmul(pset1<Packet>(Scalar(EIGEN_LN2)), x));
1860 }
1861};
1862
1863// This specialization uses a faster algorithm to compute exp2(x) for floats
1864// in [-0.5;0.5] with a relative accuracy of 1 ulp.
1865// The minimax polynomial used was calculated using the Sollya tool.
1866// See sollya.org.
1867template <>
1868struct fast_accurate_exp2<float> {
1869 template <typename Packet>
1870 EIGEN_STRONG_INLINE Packet operator()(const Packet& x) {
1871 // This function approximates exp2(x) by a degree 6 polynomial of the form
1872 // Q(x) = 1 + x * (C + x * P(x)), where the degree 4 polynomial P(x) is evaluated in
1873 // single precision, and the remaining steps are evaluated with extra precision using
1874 // double word arithmetic. C is an extra precise constant stored as a double word.
1875 //
1876 // The polynomial coefficients were calculated using Sollya commands:
1877 // > n = 6;
1878 // > f = 2^x;
1879 // > interval = [-0.5;0.5];
1880 // > p = fpminimax(f,n,[|1,double,single...|],interval,relative,floating);
1881
1882 const Packet p4 = pset1<Packet>(1.539513905e-4f);
1883 const Packet p3 = pset1<Packet>(1.340007293e-3f);
1884 const Packet p2 = pset1<Packet>(9.618283249e-3f);
1885 const Packet p1 = pset1<Packet>(5.550328270e-2f);
1886 const Packet p0 = pset1<Packet>(0.2402264923f);
1887
1888 const Packet C_hi = pset1<Packet>(0.6931471825f);
1889 const Packet C_lo = pset1<Packet>(2.36836577e-08f);
1890 const Packet one = pset1<Packet>(1.0f);
1891
1892 // Evaluate P(x) in working precision.
1893 // We evaluate even and odd parts of the polynomial separately
1894 // to gain some instruction level parallelism.
1895 Packet x2 = pmul(x, x);
1896 Packet p_even = pmadd(p4, x2, p2);
1897 Packet p_odd = pmadd(p3, x2, p1);
1898 p_even = pmadd(p_even, x2, p0);
1899 Packet p = pmadd(p_odd, x, p_even);
1900
1901 // Evaluate the remaining terms of Q(x) with extra precision using
1902 // double word arithmetic.
1903 Packet p_hi, p_lo;
1904 // x * p(x)
1905 twoprod(p, x, p_hi, p_lo);
1906 // C + x * p(x)
1907 Packet q1_hi, q1_lo;
1908 twosum(p_hi, p_lo, C_hi, C_lo, q1_hi, q1_lo);
1909 // x * (C + x * p(x))
1910 Packet q2_hi, q2_lo;
1911 twoprod(q1_hi, q1_lo, x, q2_hi, q2_lo);
1912 // 1 + x * (C + x * p(x))
1913 Packet q3_hi, q3_lo;
1914 // Since |q2_hi| <= sqrt(2)-1 < 1, we can use fast_twosum
1915 // for adding it to unity here.
1916 fast_twosum(one, q2_hi, q3_hi, q3_lo);
1917 return padd(q3_hi, padd(q2_lo, q3_lo));
1918 }
1919};
1920
1921// in [-0.5;0.5] with a relative accuracy of 1 ulp.
1922// The minimax polynomial used was calculated using the Sollya tool.
1923// See sollya.org.
1924template <>
1925struct fast_accurate_exp2<double> {
1926 template <typename Packet>
1927 EIGEN_STRONG_INLINE Packet operator()(const Packet& x) {
1928 // This function approximates exp2(x) by a degree 10 polynomial of the form
1929 // Q(x) = 1 + x * (C + x * P(x)), where the degree 8 polynomial P(x) is evaluated in
1930 // single precision, and the remaining steps are evaluated with extra precision using
1931 // double word arithmetic. C is an extra precise constant stored as a double word.
1932 //
1933 // The polynomial coefficients were calculated using Sollya commands:
1934 // > n = 11;
1935 // > f = 2^x;
1936 // > interval = [-0.5;0.5];
1937 // > p = fpminimax(f,n,[|1,DD,double...|],interval,relative,floating);
1938
1939 const Packet p9 = pset1<Packet>(4.431642109085495276e-10);
1940 const Packet p8 = pset1<Packet>(7.073829923303358410e-9);
1941 const Packet p7 = pset1<Packet>(1.017822306737031311e-7);
1942 const Packet p6 = pset1<Packet>(1.321543498017646657e-6);
1943 const Packet p5 = pset1<Packet>(1.525273342728892877e-5);
1944 const Packet p4 = pset1<Packet>(1.540353045780084423e-4);
1945 const Packet p3 = pset1<Packet>(1.333355814685869807e-3);
1946 const Packet p2 = pset1<Packet>(9.618129107593478832e-3);
1947 const Packet p1 = pset1<Packet>(5.550410866481961247e-2);
1948 const Packet p0 = pset1<Packet>(0.240226506959101332);
1949 const Packet C_hi = pset1<Packet>(0.693147180559945286);
1950 const Packet C_lo = pset1<Packet>(4.81927865669806721e-17);
1951 const Packet one = pset1<Packet>(1.0);
1952
1953 // Evaluate P(x) in working precision.
1954 // We evaluate even and odd parts of the polynomial separately
1955 // to gain some instruction level parallelism.
1956 Packet x2 = pmul(x, x);
1957 Packet p_even = pmadd(p8, x2, p6);
1958 Packet p_odd = pmadd(p9, x2, p7);
1959 p_even = pmadd(p_even, x2, p4);
1960 p_odd = pmadd(p_odd, x2, p5);
1961 p_even = pmadd(p_even, x2, p2);
1962 p_odd = pmadd(p_odd, x2, p3);
1963 p_even = pmadd(p_even, x2, p0);
1964 p_odd = pmadd(p_odd, x2, p1);
1965 Packet p = pmadd(p_odd, x, p_even);
1966
1967 // Evaluate the remaining terms of Q(x) with extra precision using
1968 // double word arithmetic.
1969 Packet p_hi, p_lo;
1970 // x * p(x)
1971 twoprod(p, x, p_hi, p_lo);
1972 // C + x * p(x)
1973 Packet q1_hi, q1_lo;
1974 twosum(p_hi, p_lo, C_hi, C_lo, q1_hi, q1_lo);
1975 // x * (C + x * p(x))
1976 Packet q2_hi, q2_lo;
1977 twoprod(q1_hi, q1_lo, x, q2_hi, q2_lo);
1978 // 1 + x * (C + x * p(x))
1979 Packet q3_hi, q3_lo;
1980 // Since |q2_hi| <= sqrt(2)-1 < 1, we can use fast_twosum
1981 // for adding it to unity here.
1982 fast_twosum(one, q2_hi, q3_hi, q3_lo);
1983 return padd(q3_hi, padd(q2_lo, q3_lo));
1984 }
1985};
1986
1987// This function implements the non-trivial case of pow(x,y) where x is
1988// positive and y is (possibly) non-integer.
1989// Formally, pow(x,y) = exp2(y * log2(x)), where exp2(x) is shorthand for 2^x.
1990// TODO(rmlarsen): We should probably add this as a packet up 'ppow', to make it
1991// easier to specialize or turn off for specific types and/or backends.x
1992template <typename Packet>
1993EIGEN_STRONG_INLINE Packet generic_pow_impl(const Packet& x, const Packet& y) {
1994 typedef typename unpacket_traits<Packet>::type Scalar;
1995 // Split x into exponent e_x and mantissa m_x.
1996 Packet e_x;
1997 Packet m_x = pfrexp(x, e_x);
1998
1999 // Adjust m_x to lie in [1/sqrt(2):sqrt(2)] to minimize absolute error in log2(m_x).
2000 EIGEN_CONSTEXPR Scalar sqrt_half = Scalar(0.70710678118654752440);
2001 const Packet m_x_scale_mask = pcmp_lt(m_x, pset1<Packet>(sqrt_half));
2002 m_x = pselect(m_x_scale_mask, pmul(pset1<Packet>(Scalar(2)), m_x), m_x);
2003 e_x = pselect(m_x_scale_mask, psub(e_x, pset1<Packet>(Scalar(1))), e_x);
2004
2005 // Compute log2(m_x) with 6 extra bits of accuracy.
2006 Packet rx_hi, rx_lo;
2007 accurate_log2<Scalar>()(m_x, rx_hi, rx_lo);
2008
2009 // Compute the two terms {y * e_x, y * r_x} in f = y * log2(x) with doubled
2010 // precision using double word arithmetic.
2011 Packet f1_hi, f1_lo, f2_hi, f2_lo;
2012 twoprod(e_x, y, f1_hi, f1_lo);
2013 twoprod(rx_hi, rx_lo, y, f2_hi, f2_lo);
2014 // Sum the two terms in f using double word arithmetic. We know
2015 // that |e_x| > |log2(m_x)|, except for the case where e_x==0.
2016 // This means that we can use fast_twosum(f1,f2).
2017 // In the case e_x == 0, e_x * y = f1 = 0, so we don't lose any
2018 // accuracy by violating the assumption of fast_twosum, because
2019 // it's a no-op.
2020 Packet f_hi, f_lo;
2021 fast_twosum(f1_hi, f1_lo, f2_hi, f2_lo, f_hi, f_lo);
2022
2023 // Split f into integer and fractional parts.
2024 Packet n_z, r_z;
2025 absolute_split(f_hi, n_z, r_z);
2026 r_z = padd(r_z, f_lo);
2027 Packet n_r;
2028 absolute_split(r_z, n_r, r_z);
2029 n_z = padd(n_z, n_r);
2030
2031 // We now have an accurate split of f = n_z + r_z and can compute
2032 // x^y = 2**{n_z + r_z) = exp2(r_z) * 2**{n_z}.
2033 // Since r_z is in [-0.5;0.5], we compute the first factor to high accuracy
2034 // using a specialized algorithm. Multiplication by the second factor can
2035 // be done exactly using pldexp(), since it is an integer power of 2.
2036 const Packet e_r = fast_accurate_exp2<Scalar>()(r_z);
2037 return pldexp(e_r, n_z);
2038}
2039
2040// Generic implementation of pow(x,y).
2041template <typename Packet>
2042EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet generic_pow(const Packet& x, const Packet& y) {
2043 typedef typename unpacket_traits<Packet>::type Scalar;
2044
2045 const Packet cst_pos_inf = pset1<Packet>(NumTraits<Scalar>::infinity());
2046 const Packet cst_neg_inf = pset1<Packet>(-NumTraits<Scalar>::infinity());
2047 const Packet cst_zero = pset1<Packet>(Scalar(0));
2048 const Packet cst_one = pset1<Packet>(Scalar(1));
2049 const Packet cst_nan = pset1<Packet>(NumTraits<Scalar>::quiet_NaN());
2050
2051 const Packet abs_x = pabs(x);
2052 // Predicates for sign and magnitude of x.
2053 const Packet abs_x_is_zero = pcmp_eq(abs_x, cst_zero);
2054 const Packet x_has_signbit = psignbit(x);
2055 const Packet x_is_neg = pandnot(x_has_signbit, abs_x_is_zero);
2056 const Packet x_is_neg_zero = pand(x_has_signbit, abs_x_is_zero);
2057 const Packet abs_x_is_inf = pcmp_eq(abs_x, cst_pos_inf);
2058 const Packet abs_x_is_one = pcmp_eq(abs_x, cst_one);
2059 const Packet abs_x_is_gt_one = pcmp_lt(cst_one, abs_x);
2060 const Packet abs_x_is_lt_one = pcmp_lt(abs_x, cst_one);
2061 const Packet x_is_one = pandnot(abs_x_is_one, x_is_neg);
2062 const Packet x_is_neg_one = pand(abs_x_is_one, x_is_neg);
2063 const Packet x_is_nan = pisnan(x);
2064
2065 // Predicates for sign and magnitude of y.
2066 const Packet abs_y = pabs(y);
2067 const Packet y_is_one = pcmp_eq(y, cst_one);
2068 const Packet abs_y_is_zero = pcmp_eq(abs_y, cst_zero);
2069 const Packet y_is_neg = pcmp_lt(y, cst_zero);
2070 const Packet y_is_pos = pandnot(ptrue(y), por(abs_y_is_zero, y_is_neg));
2071 const Packet y_is_nan = pisnan(y);
2072 const Packet abs_y_is_inf = pcmp_eq(abs_y, cst_pos_inf);
2073 EIGEN_CONSTEXPR Scalar huge_exponent =
2074 (NumTraits<Scalar>::max_exponent() * Scalar(EIGEN_LN2)) / NumTraits<Scalar>::epsilon();
2075 const Packet abs_y_is_huge = pcmp_le(pset1<Packet>(huge_exponent), pabs(y));
2076
2077 // Predicates for whether y is integer and/or even.
2078 const Packet y_is_int = pcmp_eq(pfloor(y), y);
2079 const Packet y_div_2 = pmul(y, pset1<Packet>(Scalar(0.5)));
2080 const Packet y_is_even = pcmp_eq(pround(y_div_2), y_div_2);
2081
2082 // Predicates encoding special cases for the value of pow(x,y)
2083 const Packet invalid_negative_x = pandnot(pandnot(pandnot(x_is_neg, abs_x_is_inf), y_is_int), abs_y_is_inf);
2084 const Packet pow_is_nan = por(invalid_negative_x, por(x_is_nan, y_is_nan));
2085 const Packet pow_is_one =
2086 por(por(x_is_one, abs_y_is_zero), pand(x_is_neg_one, por(abs_y_is_inf, pandnot(y_is_even, invalid_negative_x))));
2087 const Packet pow_is_zero = por(por(por(pand(abs_x_is_zero, y_is_pos), pand(abs_x_is_inf, y_is_neg)),
2088 pand(pand(abs_x_is_lt_one, abs_y_is_huge), y_is_pos)),
2089 pand(pand(abs_x_is_gt_one, abs_y_is_huge), y_is_neg));
2090 const Packet pow_is_inf = por(por(por(pand(abs_x_is_zero, y_is_neg), pand(abs_x_is_inf, y_is_pos)),
2091 pand(pand(abs_x_is_lt_one, abs_y_is_huge), y_is_neg)),
2092 pand(pand(abs_x_is_gt_one, abs_y_is_huge), y_is_pos));
2093 const Packet pow_is_neg_zero = pand(pandnot(y_is_int, y_is_even),
2094 por(pand(y_is_neg, pand(abs_x_is_inf, x_is_neg)), pand(y_is_pos, x_is_neg_zero)));
2095 const Packet inf_val =
2096 pselect(pandnot(pand(por(pand(abs_x_is_inf, x_is_neg), pand(x_is_neg_zero, y_is_neg)), y_is_int), y_is_even),
2097 cst_neg_inf, cst_pos_inf);
2098 // General computation of pow(x,y) for positive x or negative x and integer y.
2099 const Packet negate_pow_abs = pandnot(x_is_neg, y_is_even);
2100 const Packet pow_abs = generic_pow_impl(abs_x, y);
2101 return pselect(y_is_one, x,
2102 pselect(pow_is_one, cst_one,
2103 pselect(pow_is_nan, cst_nan,
2104 pselect(pow_is_inf, inf_val,
2105 pselect(pow_is_neg_zero, pnegate(cst_zero),
2106 pselect(pow_is_zero, cst_zero,
2107 pselect(negate_pow_abs, pnegate(pow_abs), pow_abs)))))));
2108}
2109
2110/* polevl (modified for Eigen)
2111 *
2112 * Evaluate polynomial
2113 *
2114 *
2115 *
2116 * SYNOPSIS:
2117 *
2118 * int N;
2119 * Scalar x, y, coef[N+1];
2120 *
2121 * y = polevl<decltype(x), N>( x, coef);
2122 *
2123 *
2124 *
2125 * DESCRIPTION:
2126 *
2127 * Evaluates polynomial of degree N:
2128 *
2129 * 2 N
2130 * y = C + C x + C x +...+ C x
2131 * 0 1 2 N
2132 *
2133 * Coefficients are stored in reverse order:
2134 *
2135 * coef[0] = C , ..., coef[N] = C .
2136 * N 0
2137 *
2138 * The function p1evl() assumes that coef[N] = 1.0 and is
2139 * omitted from the array. Its calling arguments are
2140 * otherwise the same as polevl().
2141 *
2142 *
2143 * The Eigen implementation is templatized. For best speed, store
2144 * coef as a const array (constexpr), e.g.
2145 *
2146 * const double coef[] = {1.0, 2.0, 3.0, ...};
2147 *
2148 */
2149template <typename Packet, int N>
2150struct ppolevl {
2151 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run(const Packet& x,
2152 const typename unpacket_traits<Packet>::type coeff[]) {
2153 EIGEN_STATIC_ASSERT((N > 0), YOU_MADE_A_PROGRAMMING_MISTAKE);
2154 return pmadd(ppolevl<Packet, N - 1>::run(x, coeff), x, pset1<Packet>(coeff[N]));
2155 }
2156};
2157
2158template <typename Packet>
2159struct ppolevl<Packet, 0> {
2160 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run(const Packet& x,
2161 const typename unpacket_traits<Packet>::type coeff[]) {
2162 EIGEN_UNUSED_VARIABLE(x);
2163 return pset1<Packet>(coeff[0]);
2164 }
2165};
2166
2167/* chbevl (modified for Eigen)
2168 *
2169 * Evaluate Chebyshev series
2170 *
2171 *
2172 *
2173 * SYNOPSIS:
2174 *
2175 * int N;
2176 * Scalar x, y, coef[N], chebevl();
2177 *
2178 * y = chbevl( x, coef, N );
2179 *
2180 *
2181 *
2182 * DESCRIPTION:
2183 *
2184 * Evaluates the series
2185 *
2186 * N-1
2187 * - '
2188 * y = > coef[i] T (x/2)
2189 * - i
2190 * i=0
2191 *
2192 * of Chebyshev polynomials Ti at argument x/2.
2193 *
2194 * Coefficients are stored in reverse order, i.e. the zero
2195 * order term is last in the array. Note N is the number of
2196 * coefficients, not the order.
2197 *
2198 * If coefficients are for the interval a to b, x must
2199 * have been transformed to x -> 2(2x - b - a)/(b-a) before
2200 * entering the routine. This maps x from (a, b) to (-1, 1),
2201 * over which the Chebyshev polynomials are defined.
2202 *
2203 * If the coefficients are for the inverted interval, in
2204 * which (a, b) is mapped to (1/b, 1/a), the transformation
2205 * required is x -> 2(2ab/x - b - a)/(b-a). If b is infinity,
2206 * this becomes x -> 4a/x - 1.
2207 *
2208 *
2209 *
2210 * SPEED:
2211 *
2212 * Taking advantage of the recurrence properties of the
2213 * Chebyshev polynomials, the routine requires one more
2214 * addition per loop than evaluating a nested polynomial of
2215 * the same degree.
2216 *
2217 */
2218
2219template <typename Packet, int N>
2220struct pchebevl {
2221 EIGEN_DEVICE_FUNC static EIGEN_STRONG_INLINE Packet run(Packet x,
2222 const typename unpacket_traits<Packet>::type coef[]) {
2223 typedef typename unpacket_traits<Packet>::type Scalar;
2224 Packet b0 = pset1<Packet>(coef[0]);
2225 Packet b1 = pset1<Packet>(static_cast<Scalar>(0.f));
2226 Packet b2;
2227
2228 for (int i = 1; i < N; i++) {
2229 b2 = b1;
2230 b1 = b0;
2231 b0 = psub(pmadd(x, b1, pset1<Packet>(coef[i])), b2);
2232 }
2233
2234 return pmul(pset1<Packet>(static_cast<Scalar>(0.5f)), psub(b0, b2));
2235 }
2236};
2237
2238namespace unary_pow {
2239
2240template <typename ScalarExponent, bool IsInteger = NumTraits<ScalarExponent>::IsInteger>
2241struct exponent_helper {
2242 using safe_abs_type = ScalarExponent;
2243 static constexpr ScalarExponent one_half = ScalarExponent(0.5);
2244 // these routines assume that exp is an integer stored as a floating point type
2245 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ScalarExponent safe_abs(const ScalarExponent& exp) {
2246 return numext::abs(exp);
2247 }
2248 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE bool is_odd(const ScalarExponent& exp) {
2249 eigen_assert(((numext::isfinite)(exp) && exp == numext::floor(exp)) && "exp must be an integer");
2250 ScalarExponent exp_div_2 = exp * one_half;
2251 ScalarExponent floor_exp_div_2 = numext::floor(exp_div_2);
2252 return exp_div_2 != floor_exp_div_2;
2253 }
2254 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE ScalarExponent floor_div_two(const ScalarExponent& exp) {
2255 ScalarExponent exp_div_2 = exp * one_half;
2256 return numext::floor(exp_div_2);
2257 }
2258};
2259
2260template <typename ScalarExponent>
2261struct exponent_helper<ScalarExponent, true> {
2262 // if `exp` is a signed integer type, cast it to its unsigned counterpart to safely store its absolute value
2263 // consider the (rare) case where `exp` is an int32_t: abs(-2147483648) != 2147483648
2264 using safe_abs_type = typename numext::get_integer_by_size<sizeof(ScalarExponent)>::unsigned_type;
2265 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE safe_abs_type safe_abs(const ScalarExponent& exp) {
2266 ScalarExponent mask = numext::signbit(exp);
2267 safe_abs_type result = safe_abs_type(exp ^ mask);
2268 return result + safe_abs_type(ScalarExponent(1) & mask);
2269 }
2270 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE bool is_odd(const safe_abs_type& exp) {
2271 return exp % safe_abs_type(2) != safe_abs_type(0);
2272 }
2273 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE safe_abs_type floor_div_two(const safe_abs_type& exp) {
2274 return exp >> safe_abs_type(1);
2275 }
2276};
2277
2278template <typename Packet, typename ScalarExponent,
2279 bool ReciprocateIfExponentIsNegative =
2280 !NumTraits<typename unpacket_traits<Packet>::type>::IsInteger && NumTraits<ScalarExponent>::IsSigned>
2281struct reciprocate {
2282 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run(const Packet& x, const ScalarExponent& exponent) {
2283 using Scalar = typename unpacket_traits<Packet>::type;
2284 const Packet cst_pos_one = pset1<Packet>(Scalar(1));
2285 return exponent < 0 ? pdiv(cst_pos_one, x) : x;
2286 }
2287};
2288
2289template <typename Packet, typename ScalarExponent>
2290struct reciprocate<Packet, ScalarExponent, false> {
2291 // pdiv not defined, nor necessary for integer base types
2292 // if the exponent is unsigned, then the exponent cannot be negative
2293 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run(const Packet& x, const ScalarExponent&) { return x; }
2294};
2295
2296template <typename Packet, typename ScalarExponent>
2297EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet int_pow(const Packet& x, const ScalarExponent& exponent) {
2298 using Scalar = typename unpacket_traits<Packet>::type;
2299 using ExponentHelper = exponent_helper<ScalarExponent>;
2300 using AbsExponentType = typename ExponentHelper::safe_abs_type;
2301 const Packet cst_pos_one = pset1<Packet>(Scalar(1));
2302 if (exponent == ScalarExponent(0)) return cst_pos_one;
2303
2304 Packet result = reciprocate<Packet, ScalarExponent>::run(x, exponent);
2305 Packet y = cst_pos_one;
2306 AbsExponentType m = ExponentHelper::safe_abs(exponent);
2307
2308 while (m > 1) {
2309 bool odd = ExponentHelper::is_odd(m);
2310 if (odd) y = pmul(y, result);
2311 result = pmul(result, result);
2312 m = ExponentHelper::floor_div_two(m);
2313 }
2314
2315 return pmul(y, result);
2316}
2317
2318template <typename Packet>
2319EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet gen_pow(const Packet& x,
2320 const typename unpacket_traits<Packet>::type& exponent) {
2321 const Packet exponent_packet = pset1<Packet>(exponent);
2322 return generic_pow_impl(x, exponent_packet);
2323}
2324
2325template <typename Packet, typename ScalarExponent>
2326EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet handle_nonint_nonint_errors(const Packet& x, const Packet& powx,
2327 const ScalarExponent& exponent) {
2328 using Scalar = typename unpacket_traits<Packet>::type;
2329
2330 // non-integer base and exponent case
2331
2332 const Scalar pos_zero = Scalar(0);
2333 const Scalar all_ones = ptrue<Scalar>(Scalar());
2334 const Scalar pos_one = Scalar(1);
2335 const Scalar pos_inf = NumTraits<Scalar>::infinity();
2336
2337 const Packet cst_pos_zero = pzero(x);
2338 const Packet cst_pos_one = pset1<Packet>(pos_one);
2339 const Packet cst_pos_inf = pset1<Packet>(pos_inf);
2340
2341 const bool exponent_is_not_fin = !(numext::isfinite)(exponent);
2342 const bool exponent_is_neg = exponent < ScalarExponent(0);
2343 const bool exponent_is_pos = exponent > ScalarExponent(0);
2344
2345 const Packet exp_is_not_fin = pset1<Packet>(exponent_is_not_fin ? all_ones : pos_zero);
2346 const Packet exp_is_neg = pset1<Packet>(exponent_is_neg ? all_ones : pos_zero);
2347 const Packet exp_is_pos = pset1<Packet>(exponent_is_pos ? all_ones : pos_zero);
2348 const Packet exp_is_inf = pand(exp_is_not_fin, por(exp_is_neg, exp_is_pos));
2349 const Packet exp_is_nan = pandnot(exp_is_not_fin, por(exp_is_neg, exp_is_pos));
2350
2351 const Packet x_is_le_zero = pcmp_le(x, cst_pos_zero);
2352 const Packet x_is_ge_zero = pcmp_le(cst_pos_zero, x);
2353 const Packet x_is_zero = pand(x_is_le_zero, x_is_ge_zero);
2354
2355 const Packet abs_x = pabs(x);
2356 const Packet abs_x_is_le_one = pcmp_le(abs_x, cst_pos_one);
2357 const Packet abs_x_is_ge_one = pcmp_le(cst_pos_one, abs_x);
2358 const Packet abs_x_is_inf = pcmp_eq(abs_x, cst_pos_inf);
2359 const Packet abs_x_is_one = pand(abs_x_is_le_one, abs_x_is_ge_one);
2360
2361 Packet pow_is_inf_if_exp_is_neg = por(x_is_zero, pand(abs_x_is_le_one, exp_is_inf));
2362 Packet pow_is_inf_if_exp_is_pos = por(abs_x_is_inf, pand(abs_x_is_ge_one, exp_is_inf));
2363 Packet pow_is_one = pand(abs_x_is_one, por(exp_is_inf, x_is_ge_zero));
2364
2365 Packet result = powx;
2366 result = por(x_is_le_zero, result);
2367 result = pselect(pow_is_inf_if_exp_is_neg, pand(cst_pos_inf, exp_is_neg), result);
2368 result = pselect(pow_is_inf_if_exp_is_pos, pand(cst_pos_inf, exp_is_pos), result);
2369 result = por(exp_is_nan, result);
2370 result = pselect(pow_is_one, cst_pos_one, result);
2371 return result;
2372}
2373
2374template <typename Packet, typename ScalarExponent,
2375 std::enable_if_t<NumTraits<typename unpacket_traits<Packet>::type>::IsSigned, bool> = true>
2376EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet handle_negative_exponent(const Packet& x, const ScalarExponent& exponent) {
2377 using Scalar = typename unpacket_traits<Packet>::type;
2378
2379 // singed integer base, signed integer exponent case
2380
2381 // This routine handles negative exponents.
2382 // The return value is either 0, 1, or -1.
2383
2384 const Scalar pos_zero = Scalar(0);
2385 const Scalar all_ones = ptrue<Scalar>(Scalar());
2386 const Scalar pos_one = Scalar(1);
2387
2388 const Packet cst_pos_one = pset1<Packet>(pos_one);
2389
2390 const bool exponent_is_odd = exponent % ScalarExponent(2) != ScalarExponent(0);
2391
2392 const Packet exp_is_odd = pset1<Packet>(exponent_is_odd ? all_ones : pos_zero);
2393
2394 const Packet abs_x = pabs(x);
2395 const Packet abs_x_is_one = pcmp_eq(abs_x, cst_pos_one);
2396
2397 Packet result = pselect(exp_is_odd, x, abs_x);
2398 result = pand(abs_x_is_one, result);
2399 return result;
2400}
2401
2402template <typename Packet, typename ScalarExponent,
2403 std::enable_if_t<!NumTraits<typename unpacket_traits<Packet>::type>::IsSigned, bool> = true>
2404EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet handle_negative_exponent(const Packet& x, const ScalarExponent&) {
2405 using Scalar = typename unpacket_traits<Packet>::type;
2406
2407 // unsigned integer base, signed integer exponent case
2408
2409 // This routine handles negative exponents.
2410 // The return value is either 0 or 1
2411
2412 const Scalar pos_one = Scalar(1);
2413
2414 const Packet cst_pos_one = pset1<Packet>(pos_one);
2415
2416 const Packet x_is_one = pcmp_eq(x, cst_pos_one);
2417
2418 return pand(x_is_one, x);
2419}
2420
2421} // end namespace unary_pow
2422
2423template <typename Packet, typename ScalarExponent,
2424 bool BaseIsIntegerType = NumTraits<typename unpacket_traits<Packet>::type>::IsInteger,
2425 bool ExponentIsIntegerType = NumTraits<ScalarExponent>::IsInteger,
2426 bool ExponentIsSigned = NumTraits<ScalarExponent>::IsSigned>
2427struct unary_pow_impl;
2428
2429template <typename Packet, typename ScalarExponent, bool ExponentIsSigned>
2430struct unary_pow_impl<Packet, ScalarExponent, false, false, ExponentIsSigned> {
2431 typedef typename unpacket_traits<Packet>::type Scalar;
2432 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run(const Packet& x, const ScalarExponent& exponent) {
2433 const bool exponent_is_integer = (numext::isfinite)(exponent) && numext::round(exponent) == exponent;
2434 if (exponent_is_integer) {
2435 return unary_pow::int_pow(x, exponent);
2436 } else {
2437 Packet result = unary_pow::gen_pow(x, exponent);
2438 result = unary_pow::handle_nonint_nonint_errors(x, result, exponent);
2439 return result;
2440 }
2441 }
2442};
2443
2444template <typename Packet, typename ScalarExponent, bool ExponentIsSigned>
2445struct unary_pow_impl<Packet, ScalarExponent, false, true, ExponentIsSigned> {
2446 typedef typename unpacket_traits<Packet>::type Scalar;
2447 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run(const Packet& x, const ScalarExponent& exponent) {
2448 return unary_pow::int_pow(x, exponent);
2449 }
2450};
2451
2452template <typename Packet, typename ScalarExponent>
2453struct unary_pow_impl<Packet, ScalarExponent, true, true, true> {
2454 typedef typename unpacket_traits<Packet>::type Scalar;
2455 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run(const Packet& x, const ScalarExponent& exponent) {
2456 if (exponent < ScalarExponent(0)) {
2457 return unary_pow::handle_negative_exponent(x, exponent);
2458 } else {
2459 return unary_pow::int_pow(x, exponent);
2460 }
2461 }
2462};
2463
2464template <typename Packet, typename ScalarExponent>
2465struct unary_pow_impl<Packet, ScalarExponent, true, true, false> {
2466 typedef typename unpacket_traits<Packet>::type Scalar;
2467 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run(const Packet& x, const ScalarExponent& exponent) {
2468 return unary_pow::int_pow(x, exponent);
2469 }
2470};
2471
2472template <typename Packet>
2473EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet generic_rint(const Packet& a) {
2474 using Scalar = typename unpacket_traits<Packet>::type;
2475 using IntType = typename numext::get_integer_by_size<sizeof(Scalar)>::signed_type;
2476 // Adds and subtracts signum(a) * 2^kMantissaBits to force rounding.
2477 const IntType kLimit = IntType(1) << (NumTraits<Scalar>::digits() - 1);
2478 const Packet cst_limit = pset1<Packet>(static_cast<Scalar>(kLimit));
2479 Packet abs_a = pabs(a);
2480 Packet sign_a = pandnot(a, abs_a);
2481 Packet rint_a = padd(abs_a, cst_limit);
2482 // Don't compile-away addition and subtraction.
2483 EIGEN_OPTIMIZATION_BARRIER(rint_a);
2484 rint_a = psub(rint_a, cst_limit);
2485 rint_a = por(rint_a, sign_a);
2486 // If greater than limit (or NaN), simply return a.
2487 Packet mask = pcmp_lt(abs_a, cst_limit);
2488 Packet result = pselect(mask, rint_a, a);
2489 return result;
2490}
2491
2492template <typename Packet>
2493EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet generic_floor(const Packet& a) {
2494 using Scalar = typename unpacket_traits<Packet>::type;
2495 const Packet cst_1 = pset1<Packet>(Scalar(1));
2496 Packet rint_a = generic_rint(a);
2497 // if a < rint(a), then rint(a) == ceil(a)
2498 Packet mask = pcmp_lt(a, rint_a);
2499 Packet offset = pand(cst_1, mask);
2500 Packet result = psub(rint_a, offset);
2501 return result;
2502}
2503
2504template <typename Packet>
2505EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet generic_ceil(const Packet& a) {
2506 using Scalar = typename unpacket_traits<Packet>::type;
2507 const Packet cst_1 = pset1<Packet>(Scalar(1));
2508 Packet rint_a = generic_rint(a);
2509 // if rint(a) < a, then rint(a) == floor(a)
2510 Packet mask = pcmp_lt(rint_a, a);
2511 Packet offset = pand(cst_1, mask);
2512 Packet result = padd(rint_a, offset);
2513 return result;
2514}
2515
2516template <typename Packet>
2517EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet generic_trunc(const Packet& a) {
2518 Packet abs_a = pabs(a);
2519 Packet sign_a = pandnot(a, abs_a);
2520 Packet floor_abs_a = generic_floor(abs_a);
2521 Packet result = por(floor_abs_a, sign_a);
2522 return result;
2523}
2524
2525template <typename Packet>
2526EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet generic_round(const Packet& a) {
2527 using Scalar = typename unpacket_traits<Packet>::type;
2528 const Packet cst_half = pset1<Packet>(Scalar(0.5));
2529 const Packet cst_1 = pset1<Packet>(Scalar(1));
2530 Packet abs_a = pabs(a);
2531 Packet sign_a = pandnot(a, abs_a);
2532 Packet floor_abs_a = generic_floor(abs_a);
2533 Packet diff = psub(abs_a, floor_abs_a);
2534 Packet mask = pcmp_le(cst_half, diff);
2535 Packet offset = pand(cst_1, mask);
2536 Packet result = padd(floor_abs_a, offset);
2537 result = por(result, sign_a);
2538 return result;
2539}
2540
2541template <typename Packet>
2542struct nearest_integer_packetop_impl<Packet, /*IsScalar*/ false, /*IsInteger*/ false> {
2543 using Scalar = typename unpacket_traits<Packet>::type;
2544 static_assert(packet_traits<Scalar>::HasRound, "Generic nearest integer functions are disabled for this type.");
2545 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_floor(const Packet& x) { return generic_floor(x); }
2546 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_ceil(const Packet& x) { return generic_ceil(x); }
2547 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_rint(const Packet& x) { return generic_rint(x); }
2548 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_round(const Packet& x) { return generic_round(x); }
2549 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_trunc(const Packet& x) { return generic_trunc(x); }
2550};
2551
2552template <typename Packet>
2553struct nearest_integer_packetop_impl<Packet, /*IsScalar*/ false, /*IsInteger*/ true> {
2554 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_floor(const Packet& x) { return x; }
2555 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_ceil(const Packet& x) { return x; }
2556 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_rint(const Packet& x) { return x; }
2557 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_round(const Packet& x) { return x; }
2558 static EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE Packet run_trunc(const Packet& x) { return x; }
2559};
2560
2561} // end namespace internal
2562} // end namespace Eigen
2563
2564#endif // EIGEN_ARCH_GENERIC_PACKET_MATH_FUNCTIONS_H
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_cos_op< typename Derived::Scalar >, const Derived > cos(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_expm1_op< typename Derived::Scalar >, const Derived > expm1(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_exp_op< typename Derived::Scalar >, const Derived > exp(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sign_op< typename Derived::Scalar >, const Derived > sign(const Eigen::ArrayBase< Derived > &x)