Eigen  3.3.0
arch/ZVector/MathFunctions.h
1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2007 Julien Pommier
5 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
6 //
7 // This Source Code Form is subject to the terms of the Mozilla
8 // Public License v. 2.0. If a copy of the MPL was not distributed
9 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10 
11 /* The sin, cos, exp, and log functions of this file come from
12  * Julien Pommier's sse math library: http://gruntthepeon.free.fr/ssemath/
13  */
14 
15 #ifndef EIGEN_MATH_FUNCTIONS_ALTIVEC_H
16 #define EIGEN_MATH_FUNCTIONS_ALTIVEC_H
17 
18 namespace Eigen {
19 
20 namespace internal {
21 
22 template<> EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS EIGEN_UNUSED
23 Packet2d pexp<Packet2d>(const Packet2d& _x)
24 {
25  Packet2d x = _x;
26 
27  _EIGEN_DECLARE_CONST_Packet2d(1 , 1.0);
28  _EIGEN_DECLARE_CONST_Packet2d(2 , 2.0);
29  _EIGEN_DECLARE_CONST_Packet2d(half, 0.5);
30 
31  _EIGEN_DECLARE_CONST_Packet2d(exp_hi, 709.437);
32  _EIGEN_DECLARE_CONST_Packet2d(exp_lo, -709.436139303);
33 
34  _EIGEN_DECLARE_CONST_Packet2d(cephes_LOG2EF, 1.4426950408889634073599);
35 
36  _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_p0, 1.26177193074810590878e-4);
37  _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_p1, 3.02994407707441961300e-2);
38  _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_p2, 9.99999999999999999910e-1);
39 
40  _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q0, 3.00198505138664455042e-6);
41  _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q1, 2.52448340349684104192e-3);
42  _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q2, 2.27265548208155028766e-1);
43  _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q3, 2.00000000000000000009e0);
44 
45  _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_C1, 0.693145751953125);
46  _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_C2, 1.42860682030941723212e-6);
47 
48  Packet2d tmp, fx;
49  Packet2l emm0;
50 
51  // clamp x
52  x = pmax(pmin(x, p2d_exp_hi), p2d_exp_lo);
53  /* express exp(x) as exp(g + n*log(2)) */
54  fx = pmadd(p2d_cephes_LOG2EF, x, p2d_half);
55 
56  fx = vec_floor(fx);
57 
58  tmp = pmul(fx, p2d_cephes_exp_C1);
59  Packet2d z = pmul(fx, p2d_cephes_exp_C2);
60  x = psub(x, tmp);
61  x = psub(x, z);
62 
63  Packet2d x2 = pmul(x,x);
64 
65  Packet2d px = p2d_cephes_exp_p0;
66  px = pmadd(px, x2, p2d_cephes_exp_p1);
67  px = pmadd(px, x2, p2d_cephes_exp_p2);
68  px = pmul (px, x);
69 
70  Packet2d qx = p2d_cephes_exp_q0;
71  qx = pmadd(qx, x2, p2d_cephes_exp_q1);
72  qx = pmadd(qx, x2, p2d_cephes_exp_q2);
73  qx = pmadd(qx, x2, p2d_cephes_exp_q3);
74 
75  x = pdiv(px,psub(qx,px));
76  x = pmadd(p2d_2,x,p2d_1);
77 
78  // build 2^n
79  emm0 = vec_ctsl(fx, 0);
80 
81  static const Packet2l p2l_1023 = { 1023, 1023 };
82  static const Packet2ul p2ul_52 = { 52, 52 };
83 
84  emm0 = emm0 + p2l_1023;
85  emm0 = emm0 << reinterpret_cast<Packet2l>(p2ul_52);
86 
87  // Altivec's max & min operators just drop silent NaNs. Check NaNs in
88  // inputs and return them unmodified.
89  Packet2ul isnumber_mask = reinterpret_cast<Packet2ul>(vec_cmpeq(_x, _x));
90  return vec_sel(_x, pmax(pmul(x, reinterpret_cast<Packet2d>(emm0)), _x),
91  isnumber_mask);
92 }
93 
94 template<> EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS EIGEN_UNUSED
95 Packet2d psqrt<Packet2d>(const Packet2d& x)
96 {
97  return __builtin_s390_vfsqdb(x);
98 }
99 
100 template<> EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS EIGEN_UNUSED
101 Packet2d prsqrt<Packet2d>(const Packet2d& x) {
102  // Unfortunately we can't use the much faster mm_rqsrt_pd since it only provides an approximation.
103  return pset1<Packet2d>(1.0) / psqrt<Packet2d>(x);
104 }
105 
106 } // end namespace internal
107 
108 } // end namespace Eigen
109 
110 #endif // EIGEN_MATH_FUNCTIONS_ALTIVEC_H
Namespace containing all symbols from the Eigen library.
Definition: Core:287
Definition: Eigen_Colamd.h:50