![]() |
Eigen
3.3.3
|
00001 // This file is part of Eigen, a lightweight C++ template library 00002 // for linear algebra. 00003 // 00004 // Copyright (C) 2007 Julien Pommier 00005 // Copyright (C) 2009 Gael Guennebaud <gael.guennebaud@inria.fr> 00006 // Copyright (C) 2016 Konstantinos Margaritis <markos@freevec.org> 00007 // 00008 // This Source Code Form is subject to the terms of the Mozilla 00009 // Public License v. 2.0. If a copy of the MPL was not distributed 00010 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/. 00011 00012 /* The sin, cos, exp, and log functions of this file come from 00013 * Julien Pommier's sse math library: http://gruntthepeon.free.fr/ssemath/ 00014 */ 00015 00016 #ifndef EIGEN_MATH_FUNCTIONS_ALTIVEC_H 00017 #define EIGEN_MATH_FUNCTIONS_ALTIVEC_H 00018 00019 namespace Eigen { 00020 00021 namespace internal { 00022 00023 static _EIGEN_DECLARE_CONST_Packet2d(1 , 1.0); 00024 static _EIGEN_DECLARE_CONST_Packet2d(2 , 2.0); 00025 static _EIGEN_DECLARE_CONST_Packet2d(half, 0.5); 00026 00027 static _EIGEN_DECLARE_CONST_Packet2d(exp_hi, 709.437); 00028 static _EIGEN_DECLARE_CONST_Packet2d(exp_lo, -709.436139303); 00029 00030 static _EIGEN_DECLARE_CONST_Packet2d(cephes_LOG2EF, 1.4426950408889634073599); 00031 00032 static _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_p0, 1.26177193074810590878e-4); 00033 static _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_p1, 3.02994407707441961300e-2); 00034 static _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_p2, 9.99999999999999999910e-1); 00035 00036 static _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q0, 3.00198505138664455042e-6); 00037 static _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q1, 2.52448340349684104192e-3); 00038 static _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q2, 2.27265548208155028766e-1); 00039 static _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_q3, 2.00000000000000000009e0); 00040 00041 static _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_C1, 0.693145751953125); 00042 static _EIGEN_DECLARE_CONST_Packet2d(cephes_exp_C2, 1.42860682030941723212e-6); 00043 00044 template<> EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS EIGEN_UNUSED 00045 Packet2d pexp<Packet2d>(const Packet2d& _x) 00046 { 00047 Packet2d x = _x; 00048 00049 Packet2d tmp, fx; 00050 Packet2l emm0; 00051 00052 // clamp x 00053 x = pmax(pmin(x, p2d_exp_hi), p2d_exp_lo); 00054 /* express exp(x) as exp(g + n*log(2)) */ 00055 fx = pmadd(p2d_cephes_LOG2EF, x, p2d_half); 00056 00057 fx = vec_floor(fx); 00058 00059 tmp = pmul(fx, p2d_cephes_exp_C1); 00060 Packet2d z = pmul(fx, p2d_cephes_exp_C2); 00061 x = psub(x, tmp); 00062 x = psub(x, z); 00063 00064 Packet2d x2 = pmul(x,x); 00065 00066 Packet2d px = p2d_cephes_exp_p0; 00067 px = pmadd(px, x2, p2d_cephes_exp_p1); 00068 px = pmadd(px, x2, p2d_cephes_exp_p2); 00069 px = pmul (px, x); 00070 00071 Packet2d qx = p2d_cephes_exp_q0; 00072 qx = pmadd(qx, x2, p2d_cephes_exp_q1); 00073 qx = pmadd(qx, x2, p2d_cephes_exp_q2); 00074 qx = pmadd(qx, x2, p2d_cephes_exp_q3); 00075 00076 x = pdiv(px,psub(qx,px)); 00077 x = pmadd(p2d_2,x,p2d_1); 00078 00079 // build 2^n 00080 emm0 = vec_ctsl(fx, 0); 00081 00082 static const Packet2l p2l_1023 = { 1023, 1023 }; 00083 static const Packet2ul p2ul_52 = { 52, 52 }; 00084 00085 emm0 = emm0 + p2l_1023; 00086 emm0 = emm0 << reinterpret_cast<Packet2l>(p2ul_52); 00087 00088 // Altivec's max & min operators just drop silent NaNs. Check NaNs in 00089 // inputs and return them unmodified. 00090 Packet2ul isnumber_mask = reinterpret_cast<Packet2ul>(vec_cmpeq(_x, _x)); 00091 return vec_sel(_x, pmax(pmul(x, reinterpret_cast<Packet2d>(emm0)), _x), 00092 isnumber_mask); 00093 } 00094 00095 template<> EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS EIGEN_UNUSED 00096 Packet4f pexp<Packet4f>(const Packet4f& x) 00097 { 00098 Packet4f res; 00099 res.v4f[0] = pexp<Packet2d>(x.v4f[0]); 00100 res.v4f[1] = pexp<Packet2d>(x.v4f[1]); 00101 return res; 00102 } 00103 00104 template<> EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS EIGEN_UNUSED 00105 Packet2d psqrt<Packet2d>(const Packet2d& x) 00106 { 00107 return __builtin_s390_vfsqdb(x); 00108 } 00109 00110 template<> EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS EIGEN_UNUSED 00111 Packet4f psqrt<Packet4f>(const Packet4f& x) 00112 { 00113 Packet4f res; 00114 res.v4f[0] = psqrt<Packet2d>(x.v4f[0]); 00115 res.v4f[1] = psqrt<Packet2d>(x.v4f[1]); 00116 return res; 00117 } 00118 00119 template<> EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS EIGEN_UNUSED 00120 Packet2d prsqrt<Packet2d>(const Packet2d& x) { 00121 // Unfortunately we can't use the much faster mm_rqsrt_pd since it only provides an approximation. 00122 return pset1<Packet2d>(1.0) / psqrt<Packet2d>(x); 00123 } 00124 00125 template<> EIGEN_DEFINE_FUNCTION_ALLOWING_MULTIPLE_DEFINITIONS EIGEN_UNUSED 00126 Packet4f prsqrt<Packet4f>(const Packet4f& x) { 00127 Packet4f res; 00128 res.v4f[0] = prsqrt<Packet2d>(x.v4f[0]); 00129 res.v4f[1] = prsqrt<Packet2d>(x.v4f[1]); 00130 return res; 00131 } 00132 00133 } // end namespace internal 00134 00135 } // end namespace Eigen 00136 00137 #endif // EIGEN_MATH_FUNCTIONS_ALTIVEC_H