1
0
mirror of https://github.com/paboyle/Grid.git synced 2025-04-04 19:25:56 +01:00

Accelerator flaggina dn thrust complex for NVCC

This commit is contained in:
paboyle 2018-01-24 13:50:41 +00:00
parent 725f03e2e2
commit 0626c1e39e

View File

@ -38,142 +38,142 @@ NAMESPACE_BEGIN(Grid);
template <class scalar>
struct SqrtRealFunctor {
scalar operator()(const scalar &a) const { return sqrt(real(a)); }
accelerator scalar operator()(const scalar &a) const { return sqrt(real(a)); }
};
template <class scalar>
struct RSqrtRealFunctor {
scalar operator()(const scalar &a) const {
accelerator scalar operator()(const scalar &a) const {
return scalar(1.0 / sqrt(real(a)));
}
};
template <class scalar>
struct CosRealFunctor {
scalar operator()(const scalar &a) const { return cos(real(a)); }
accelerator scalar operator()(const scalar &a) const { return cos(real(a)); }
};
template <class scalar>
struct SinRealFunctor {
scalar operator()(const scalar &a) const { return sin(real(a)); }
accelerator scalar operator()(const scalar &a) const { return sin(real(a)); }
};
template <class scalar>
struct AcosRealFunctor {
scalar operator()(const scalar &a) const { return acos(real(a)); }
accelerator scalar operator()(const scalar &a) const { return acos(real(a)); }
};
template <class scalar>
struct AsinRealFunctor {
scalar operator()(const scalar &a) const { return asin(real(a)); }
accelerator scalar operator()(const scalar &a) const { return asin(real(a)); }
};
template <class scalar>
struct LogRealFunctor {
scalar operator()(const scalar &a) const { return log(real(a)); }
accelerator scalar operator()(const scalar &a) const { return log(real(a)); }
};
template <class scalar>
struct ExpFunctor {
scalar operator()(const scalar &a) const { return exp(a); }
accelerator scalar operator()(const scalar &a) const { return exp(a); }
};
template <class scalar>
struct NotFunctor {
scalar operator()(const scalar &a) const { return (!a); }
accelerator scalar operator()(const scalar &a) const { return (!a); }
};
template <class scalar>
struct AbsRealFunctor {
scalar operator()(const scalar &a) const { return std::abs(real(a)); }
accelerator scalar operator()(const scalar &a) const { return std::abs(real(a)); }
};
template <class scalar>
struct PowRealFunctor {
double y;
PowRealFunctor(double _y) : y(_y){};
scalar operator()(const scalar &a) const { return pow(real(a), y); }
accelerator PowRealFunctor(double _y) : y(_y){};
accelerator scalar operator()(const scalar &a) const { return pow(real(a), y); }
};
template <class scalar>
struct ModIntFunctor {
Integer y;
ModIntFunctor(Integer _y) : y(_y){};
scalar operator()(const scalar &a) const { return Integer(a) % y; }
accelerator ModIntFunctor(Integer _y) : y(_y){};
accelerator scalar operator()(const scalar &a) const { return Integer(a) % y; }
};
template <class scalar>
struct DivIntFunctor {
Integer y;
DivIntFunctor(Integer _y) : y(_y){};
scalar operator()(const scalar &a) const { return Integer(a) / y; }
accelerator DivIntFunctor(Integer _y) : y(_y){};
accelerator scalar operator()(const scalar &a) const { return Integer(a) / y; }
};
template <class scalar>
struct RealFunctor {
scalar operator()(const scalar &a) const { return std::real(a); }
accelerator scalar operator()(const scalar &a) const { return real(a); }
};
template <class scalar>
struct ImagFunctor {
scalar operator()(const scalar &a) const { return std::imag(a); }
accelerator scalar operator()(const scalar &a) const { return imag(a); }
};
template <class S, class V>
inline Grid_simd<S, V> real(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> real(const Grid_simd<S, V> &r) {
return SimdApply(RealFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> imag(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> imag(const Grid_simd<S, V> &r) {
return SimdApply(ImagFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> sqrt(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> sqrt(const Grid_simd<S, V> &r) {
return SimdApply(SqrtRealFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> rsqrt(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> rsqrt(const Grid_simd<S, V> &r) {
return SimdApply(RSqrtRealFunctor<S>(), r);
}
template <class Scalar>
inline Scalar rsqrt(const Scalar &r) {
accelerator_inline Scalar rsqrt(const Scalar &r) {
return (RSqrtRealFunctor<Scalar>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> cos(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> cos(const Grid_simd<S, V> &r) {
return SimdApply(CosRealFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> sin(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> sin(const Grid_simd<S, V> &r) {
return SimdApply(SinRealFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> acos(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> acos(const Grid_simd<S, V> &r) {
return SimdApply(AcosRealFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> asin(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> asin(const Grid_simd<S, V> &r) {
return SimdApply(AsinRealFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> log(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> log(const Grid_simd<S, V> &r) {
return SimdApply(LogRealFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> abs(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> abs(const Grid_simd<S, V> &r) {
return SimdApply(AbsRealFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> exp(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> exp(const Grid_simd<S, V> &r) {
return SimdApply(ExpFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> Not(const Grid_simd<S, V> &r) {
accelerator_inline Grid_simd<S, V> Not(const Grid_simd<S, V> &r) {
return SimdApply(NotFunctor<S>(), r);
}
template <class S, class V>
inline Grid_simd<S, V> pow(const Grid_simd<S, V> &r, double y) {
accelerator_inline Grid_simd<S, V> pow(const Grid_simd<S, V> &r, double y) {
return SimdApply(PowRealFunctor<S>(y), r);
}
template <class S, class V>
inline Grid_simd<S, V> mod(const Grid_simd<S, V> &r, Integer y) {
accelerator_inline Grid_simd<S, V> mod(const Grid_simd<S, V> &r, Integer y) {
return SimdApply(ModIntFunctor<S>(y), r);
}
template <class S, class V>
inline Grid_simd<S, V> div(const Grid_simd<S, V> &r, Integer y) {
accelerator_inline Grid_simd<S, V> div(const Grid_simd<S, V> &r, Integer y) {
return SimdApply(DivIntFunctor<S>(y), r);
}
////////////////////////////////////////////////////////////////////////////
@ -181,41 +181,41 @@ inline Grid_simd<S, V> div(const Grid_simd<S, V> &r, Integer y) {
////////////////////////////////////////////////////////////////////////////
template <class scalar>
struct AndFunctor {
scalar operator()(const scalar &x, const scalar &y) const { return x & y; }
accelerator scalar operator()(const scalar &x, const scalar &y) const { return x & y; }
};
template <class scalar>
struct OrFunctor {
scalar operator()(const scalar &x, const scalar &y) const { return x | y; }
accelerator scalar operator()(const scalar &x, const scalar &y) const { return x | y; }
};
template <class scalar>
struct AndAndFunctor {
scalar operator()(const scalar &x, const scalar &y) const { return x && y; }
accelerator scalar operator()(const scalar &x, const scalar &y) const { return x && y; }
};
template <class scalar>
struct OrOrFunctor {
scalar operator()(const scalar &x, const scalar &y) const { return x || y; }
accelerator scalar operator()(const scalar &x, const scalar &y) const { return x || y; }
};
////////////////////////////////
// Calls to simd binop functors
////////////////////////////////
template <class S, class V>
inline Grid_simd<S, V> operator&(const Grid_simd<S, V> &x,
accelerator_inline Grid_simd<S, V> operator&(const Grid_simd<S, V> &x,
const Grid_simd<S, V> &y) {
return SimdApplyBinop(AndFunctor<S>(), x, y);
}
template <class S, class V>
inline Grid_simd<S, V> operator&&(const Grid_simd<S, V> &x,
accelerator_inline Grid_simd<S, V> operator&&(const Grid_simd<S, V> &x,
const Grid_simd<S, V> &y) {
return SimdApplyBinop(AndAndFunctor<S>(), x, y);
}
template <class S, class V>
inline Grid_simd<S, V> operator|(const Grid_simd<S, V> &x,
accelerator_inline Grid_simd<S, V> operator|(const Grid_simd<S, V> &x,
const Grid_simd<S, V> &y) {
return SimdApplyBinop(OrFunctor<S>(), x, y);
}
template <class S, class V>
inline Grid_simd<S, V> operator||(const Grid_simd<S, V> &x,
accelerator_inline Grid_simd<S, V> operator||(const Grid_simd<S, V> &x,
const Grid_simd<S, V> &y) {
return SimdApplyBinop(OrOrFunctor<S>(), x, y);
}