arch/ZVector/MathFunctions.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) 2009 Gael Guennebaud <gael.guennebaud@inria.fr>
6 // Copyright (C) 2016 Konstantinos Margaritis <markos@freevec.org>
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 sin, cos, exp, and log functions of this file come from
13  * Julien Pommier's sse math library: http://gruntthepeon.free.fr/ssemath/
14  */
15 
16 #ifndef EIGEN_MATH_FUNCTIONS_ZVECTOR_H
17 #define EIGEN_MATH_FUNCTIONS_ZVECTOR_H
18 
19 #include "../../InternalHeaderCheck.h"
20 
21 namespace Eigen {
22 
23 namespace internal {
24 
25 #if !defined(__ARCH__) || (defined(__ARCH__) && __ARCH__ >= 12)
28 static EIGEN_DECLARE_CONST_Packet4i(0x7f, 0x7f);
30 
31 static EIGEN_DECLARE_CONST_Packet4f_FROM_INT(inv_mant_mask, ~0x7f800000);
32 
33 /* the smallest non denormalized float number */
34 static EIGEN_DECLARE_CONST_Packet4f_FROM_INT(min_norm_pos, 0x00800000);
35 static EIGEN_DECLARE_CONST_Packet4f_FROM_INT(minus_inf, 0xff800000); // -1.f/0.f
36 static EIGEN_DECLARE_CONST_Packet4f_FROM_INT(minus_nan, 0xffffffff);
37 
38 /* natural logarithm computed for 4 simultaneous float
39  return NaN for x <= 0
40 */
41 static EIGEN_DECLARE_CONST_Packet4f(cephes_SQRTHF, 0.707106781186547524f);
42 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_p0, 7.0376836292E-2f);
43 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_p1, - 1.1514610310E-1f);
44 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_p2, 1.1676998740E-1f);
45 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_p3, - 1.2420140846E-1f);
46 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_p4, + 1.4249322787E-1f);
47 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_p5, - 1.6668057665E-1f);
48 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_p6, + 2.0000714765E-1f);
49 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_p7, - 2.4999993993E-1f);
50 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_p8, + 3.3333331174E-1f);
51 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_q1, -2.12194440e-4f);
52 static EIGEN_DECLARE_CONST_Packet4f(cephes_log_q2, 0.693359375f);
53 
54 static EIGEN_DECLARE_CONST_Packet4f(exp_hi, 88.3762626647950f);
55 static EIGEN_DECLARE_CONST_Packet4f(exp_lo, -88.3762626647949f);
56 
57 static EIGEN_DECLARE_CONST_Packet4f(cephes_LOG2EF, 1.44269504088896341f);
58 static EIGEN_DECLARE_CONST_Packet4f(cephes_exp_C1, 0.693359375f);
59 static EIGEN_DECLARE_CONST_Packet4f(cephes_exp_C2, -2.12194440e-4f);
60 
61 static EIGEN_DECLARE_CONST_Packet4f(cephes_exp_p0, 1.9875691500E-4f);
62 static EIGEN_DECLARE_CONST_Packet4f(cephes_exp_p1, 1.3981999507E-3f);
63 static EIGEN_DECLARE_CONST_Packet4f(cephes_exp_p2, 8.3334519073E-3f);
64 static EIGEN_DECLARE_CONST_Packet4f(cephes_exp_p3, 4.1665795894E-2f);
65 static EIGEN_DECLARE_CONST_Packet4f(cephes_exp_p4, 1.6666665459E-1f);
66 static EIGEN_DECLARE_CONST_Packet4f(cephes_exp_p5, 5.0000001201E-1f);
67 #endif
68 
72 
73 static EIGEN_DECLARE_CONST_Packet2d(exp_hi, 709.437);
74 static EIGEN_DECLARE_CONST_Packet2d(exp_lo, -709.436139303);
75 
76 static EIGEN_DECLARE_CONST_Packet2d(cephes_LOG2EF, 1.4426950408889634073599);
77 
78 static EIGEN_DECLARE_CONST_Packet2d(cephes_exp_p0, 1.26177193074810590878e-4);
79 static EIGEN_DECLARE_CONST_Packet2d(cephes_exp_p1, 3.02994407707441961300e-2);
80 static EIGEN_DECLARE_CONST_Packet2d(cephes_exp_p2, 9.99999999999999999910e-1);
81 
82 static EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q0, 3.00198505138664455042e-6);
83 static EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q1, 2.52448340349684104192e-3);
84 static EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q2, 2.27265548208155028766e-1);
85 static EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q3, 2.00000000000000000009e0);
86 
87 static EIGEN_DECLARE_CONST_Packet2d(cephes_exp_C1, 0.693145751953125);
88 static EIGEN_DECLARE_CONST_Packet2d(cephes_exp_C2, 1.42860682030941723212e-6);
89 
92 {
93  Packet2d x = _x;
94 
95  Packet2d tmp, fx;
96  Packet2l emm0;
97 
98  // clamp x
99  x = pmax(pmin(x, p2d_exp_hi), p2d_exp_lo);
100  /* express exp(x) as exp(g + n*log(2)) */
101  fx = pmadd(p2d_cephes_LOG2EF, x, p2d_half);
102 
103  fx = vec_floor(fx);
104 
105  tmp = pmul(fx, p2d_cephes_exp_C1);
106  Packet2d z = pmul(fx, p2d_cephes_exp_C2);
107  x = psub(x, tmp);
108  x = psub(x, z);
109 
110  Packet2d x2 = pmul(x,x);
111 
112  Packet2d px = p2d_cephes_exp_p0;
113  px = pmadd(px, x2, p2d_cephes_exp_p1);
114  px = pmadd(px, x2, p2d_cephes_exp_p2);
115  px = pmul (px, x);
116 
117  Packet2d qx = p2d_cephes_exp_q0;
118  qx = pmadd(qx, x2, p2d_cephes_exp_q1);
119  qx = pmadd(qx, x2, p2d_cephes_exp_q2);
120  qx = pmadd(qx, x2, p2d_cephes_exp_q3);
121 
122  x = pdiv(px,psub(qx,px));
123  x = pmadd(p2d_2,x,p2d_1);
124 
125  // build 2^n
126  emm0 = vec_ctsl(fx, 0);
127 
128  static const Packet2l p2l_1023 = { 1023, 1023 };
129  static const Packet2ul p2ul_52 = { 52, 52 };
130 
131  emm0 = emm0 + p2l_1023;
132  emm0 = emm0 << reinterpret_cast<Packet2l>(p2ul_52);
133 
134  // Altivec's max & min operators just drop silent NaNs. Check NaNs in
135  // inputs and return them unmodified.
136  Packet2ul isnumber_mask = reinterpret_cast<Packet2ul>(vec_cmpeq(_x, _x));
137  return vec_sel(_x, pmax(pmul(x, reinterpret_cast<Packet2d>(emm0)), _x),
138  isnumber_mask);
139 }
140 
143 {
144 #if !defined(__ARCH__) || (defined(__ARCH__) && __ARCH__ >= 12)
145  Packet4f x = _x;
146 
147  Packet4f tmp, fx;
148  Packet4i emm0;
149 
150  // clamp x
151  x = pmax(pmin(x, p4f_exp_hi), p4f_exp_lo);
152 
153  // express exp(x) as exp(g + n*log(2))
154  fx = pmadd(x, p4f_cephes_LOG2EF, p4f_half);
155 
156  fx = pfloor(fx);
157 
158  tmp = pmul(fx, p4f_cephes_exp_C1);
159  Packet4f z = pmul(fx, p4f_cephes_exp_C2);
160  x = psub(x, tmp);
161  x = psub(x, z);
162 
163  z = pmul(x,x);
164 
165  Packet4f y = p4f_cephes_exp_p0;
166  y = pmadd(y, x, p4f_cephes_exp_p1);
167  y = pmadd(y, x, p4f_cephes_exp_p2);
168  y = pmadd(y, x, p4f_cephes_exp_p3);
169  y = pmadd(y, x, p4f_cephes_exp_p4);
170  y = pmadd(y, x, p4f_cephes_exp_p5);
171  y = pmadd(y, z, x);
172  y = padd(y, p4f_1);
173 
174  // build 2^n
175  emm0 = (Packet4i){ (int)fx[0], (int)fx[1], (int)fx[2], (int)fx[3] };
176  emm0 = emm0 + p4i_0x7f;
177  emm0 = emm0 << reinterpret_cast<Packet4i>(p4i_23);
178 
179  return pmax(pmul(y, reinterpret_cast<Packet4f>(emm0)), _x);
180 #else
181  Packet4f res;
182  res.v4f[0] = pexp<Packet2d>(_x.v4f[0]);
183  res.v4f[1] = pexp<Packet2d>(_x.v4f[1]);
184  return res;
185 #endif
186 }
187 
190 {
191  return vec_sqrt(x);
192 }
193 
196 {
197  Packet4f res;
198 #if !defined(__ARCH__) || (defined(__ARCH__) && __ARCH__ >= 12)
199  res = vec_sqrt(x);
200 #else
201  res.v4f[0] = psqrt<Packet2d>(x.v4f[0]);
202  res.v4f[1] = psqrt<Packet2d>(x.v4f[1]);
203 #endif
204  return res;
205 }
206 
209  return pset1<Packet2d>(1.0) / psqrt<Packet2d>(x);
210 }
211 
214  Packet4f res;
215 #if !defined(__ARCH__) || (defined(__ARCH__) && __ARCH__ >= 12)
217 #else
218  res.v4f[0] = prsqrt<Packet2d>(x.v4f[0]);
219  res.v4f[1] = prsqrt<Packet2d>(x.v4f[1]);
220 #endif
221  return res;
222 }
223 
224 // Hyperbolic Tangent function.
225 template <>
227 ptanh<Packet4f>(const Packet4f& x) {
229 }
230 
231 } // end namespace internal
232 
233 } // end namespace Eigen
234 
235 #endif // EIGEN_MATH_FUNCTIONS_ZVECTOR_H
#define EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS
Definition: Macros.h:892
cout<< "Here is the matrix m:"<< endl<< m<< endl;Matrix< ptrdiff_t, 3, 1 > res
Packet pmin(const Packet &a, const Packet &b)
Packet padd(const Packet &a, const Packet &b)
__vector int Packet4i
EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet2d pexp< Packet2d >(const Packet2d &_x)
EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet4f psqrt< Packet4f >(const Packet4f &x)
static EIGEN_DECLARE_CONST_Packet4f(1, 1.0f)
const Scalar & y
EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet2d prsqrt< Packet2d >(const Packet2d &x)
Packet pmax(const Packet &a, const Packet &b)
static EIGEN_DECLARE_CONST_Packet4f_FROM_INT(inv_mant_mask, ~0x7f800000)
static EIGEN_DECLARE_CONST_Packet4i(0x7f, 0x7f)
Packet4f pmadd(const Packet4f &a, const Packet4f &b, const Packet4f &c)
Packet pmul(const Packet &a, const Packet &b)
EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet4f ptanh< Packet4f >(const Packet4f &_x)
Packet psub(const Packet &a, const Packet &b)
Packet2d pset1< Packet2d >(const double &from)
T generic_fast_tanh_float(const T &a_x)
EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet2d psqrt< Packet2d >(const Packet2d &x)
Packet pdiv(const Packet &a, const Packet &b)
EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet4f pexp< Packet4f >(const Packet4f &_x)
static EIGEN_DECLARE_CONST_Packet2d(1, 1.0)
Packet4f pset1< Packet4f >(const float &from)
EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet4f prsqrt< Packet4f >(const Packet4f &x)
__vector float Packet4f
EIGEN_DECLARE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS Packet pfloor(const Packet &a)
: InteropHeaders
Definition: Core:139