11#ifndef EIGEN_MATHFUNCTIONSIMPL_H
12#define EIGEN_MATHFUNCTIONSIMPL_H
15#include "./InternalHeaderCheck.h"
35template <
typename Packet,
int Steps>
36struct generic_reciprocal_newton_step {
37 static_assert(Steps > 0,
"Steps must be at least 1.");
38 EIGEN_DEVICE_FUNC
static EIGEN_STRONG_INLINE Packet run(
const Packet& a,
const Packet& approx_a_recip) {
39 using Scalar =
typename unpacket_traits<Packet>::type;
40 const Packet two = pset1<Packet>(Scalar(2));
43 const Packet x = generic_reciprocal_newton_step<Packet, Steps - 1>::run(a, approx_a_recip);
44 const Packet tmp = pnmadd(a, x, two);
47 const Packet is_not_nan = pcmp_eq(tmp, tmp);
48 return pselect(is_not_nan, pmul(x, tmp), x);
52template <
typename Packet>
53struct generic_reciprocal_newton_step<Packet, 0> {
54 EIGEN_DEVICE_FUNC
static EIGEN_STRONG_INLINE Packet run(
const Packet& ,
const Packet& approx_rsqrt) {
74template <
typename Packet,
int Steps>
75struct generic_rsqrt_newton_step {
76 static_assert(Steps > 0,
"Steps must be at least 1.");
77 using Scalar =
typename unpacket_traits<Packet>::type;
78 EIGEN_DEVICE_FUNC
static EIGEN_STRONG_INLINE Packet run(
const Packet& a,
const Packet& approx_rsqrt) {
79 constexpr Scalar kMinusHalf = Scalar(-1) / Scalar(2);
80 const Packet cst_minus_half = pset1<Packet>(kMinusHalf);
81 const Packet cst_minus_one = pset1<Packet>(Scalar(-1));
83 Packet inv_sqrt = approx_rsqrt;
84 for (
int step = 0; step < Steps; ++step) {
88 Packet r2 = pmul(a, inv_sqrt);
89 Packet half_r = pmul(inv_sqrt, cst_minus_half);
90 Packet h_n = pmadd(r2, inv_sqrt, cst_minus_one);
91 inv_sqrt = pmadd(half_r, h_n, inv_sqrt);
98 return pselect(pisnan(inv_sqrt), approx_rsqrt, inv_sqrt);
102template <
typename Packet>
103struct generic_rsqrt_newton_step<Packet, 0> {
104 EIGEN_DEVICE_FUNC
static EIGEN_STRONG_INLINE Packet run(
const Packet& ,
const Packet& approx_rsqrt) {
124template <
typename Packet,
int Steps = 1>
125struct generic_sqrt_newton_step {
126 static_assert(Steps > 0,
"Steps must be at least 1.");
128 EIGEN_DEVICE_FUNC
static EIGEN_STRONG_INLINE Packet run(
const Packet& a,
const Packet& approx_rsqrt) {
129 using Scalar =
typename unpacket_traits<Packet>::type;
130 const Packet one_point_five = pset1<Packet>(Scalar(1.5));
131 const Packet minus_half = pset1<Packet>(Scalar(-0.5));
133 const Packet inf_mask = pcmp_eq(a, pset1<Packet>(NumTraits<Scalar>::infinity()));
134 const Packet return_a = por(pcmp_eq(a, pzero(a)), inf_mask);
138 Packet
rsqrt = pmul(approx_rsqrt, pmadd(pmul(minus_half, approx_rsqrt), pmul(a, approx_rsqrt), one_point_five));
139 for (
int step = 1; step < Steps; ++step) {
145 return pselect(return_a, a, pmul(a,
rsqrt));
149template <
typename RealScalar>
150EIGEN_DEVICE_FUNC EIGEN_STRONG_INLINE RealScalar positive_real_hypot(
const RealScalar& x,
const RealScalar& y) {
152 if ((numext::isinf)(x) || (numext::isinf)(y))
return NumTraits<RealScalar>::infinity();
153 if ((numext::isnan)(x) || (numext::isnan)(y))
return NumTraits<RealScalar>::quiet_NaN();
155 EIGEN_USING_STD(
sqrt);
157 p = numext::maxi(x, y);
158 if (numext::is_exactly_zero(p))
return RealScalar(0);
159 qp = numext::mini(y, x) / p;
160 return p *
sqrt(RealScalar(1) + qp * qp);
163template <
typename Scalar>
165 typedef typename NumTraits<Scalar>::Real RealScalar;
166 static EIGEN_DEVICE_FUNC
inline RealScalar run(
const Scalar& x,
const Scalar& y) {
167 EIGEN_USING_STD(
abs);
168 return positive_real_hypot<RealScalar>(
abs(x),
abs(y));
175EIGEN_DEVICE_FUNC std::complex<T> complex_sqrt(
const std::complex<T>& z) {
198 const T x = numext::real(z);
199 const T y = numext::imag(z);
201 const T w = numext::sqrt(T(0.5) * (numext::abs(x) + numext::hypot(x, y)));
203 return (numext::isinf)(y) ? std::complex<T>(NumTraits<T>::infinity(), y)
204 : numext::is_exactly_zero(x) ? std::complex<T>(w, y < zero ? -w : w)
205 : x > zero ? std::complex<T>(w, y / (2 * w))
206 : std::complex<T>(numext::
abs(y) / (2 * w), y < zero ? -w : w);
211EIGEN_DEVICE_FUNC std::complex<T> complex_rsqrt(
const std::complex<T>& z) {
234 const T x = numext::real(z);
235 const T y = numext::imag(z);
238 const T abs_z = numext::hypot(x, y);
239 const T w = numext::sqrt(T(0.5) * (numext::abs(x) + abs_z));
240 const T woz = w / abs_z;
242 return numext::is_exactly_zero(abs_z) ? std::complex<T>(NumTraits<T>::infinity(), NumTraits<T>::quiet_NaN())
243 : ((numext::
isinf)(x) || (numext::
isinf)(y)) ? std::complex<T>(zero, zero)
244 : numext::is_exactly_zero(x) ? std::complex<T>(woz, y < zero ? woz : -woz)
245 : x > zero ? std::complex<T>(woz, -y / (2 * w * abs_z))
246 : std::complex<T>(numext::
abs(y) / (2 * w * abs_z), y < zero ? woz : -woz);
250EIGEN_DEVICE_FUNC std::complex<T> complex_log(
const std::complex<T>& z) {
252 T a = numext::abs(z);
253 EIGEN_USING_STD(atan2);
254 T b = atan2(z.imag(), z.real());
255 return std::complex<T>(numext::log(a), b);
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_rsqrt_op< typename Derived::Scalar >, const Derived > rsqrt(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_abs_op< typename Derived::Scalar >, const Derived > abs(const Eigen::ArrayBase< Derived > &x)
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_isinf_op< typename Derived::Scalar >, const Derived > isinf(const Eigen::ArrayBase< Derived > &x)