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