mirror of
https://github.com/paboyle/Grid.git
synced 2025-06-12 20:27:06 +01:00
Conjugate residual algorithm; some more unary functions
This commit is contained in:
@ -314,9 +314,9 @@ namespace Optimization {
|
||||
template<>
|
||||
inline Grid::ComplexF Reduce<Grid::ComplexF, __m256>::operator()(__m256 in){
|
||||
__m256 v1,v2;
|
||||
Optimization::permute(v1,in,0); // sse 128; paired complex single
|
||||
Optimization::permute(v1,in,0); // avx 256; quad complex single
|
||||
v1 = _mm256_add_ps(v1,in);
|
||||
Optimization::permute(v2,v1,1); // avx 256; quad complex single
|
||||
Optimization::permute(v2,v1,1);
|
||||
v1 = _mm256_add_ps(v1,v2);
|
||||
u256f conv; conv.v = v1;
|
||||
return Grid::ComplexF(conv.f[0],conv.f[1]);
|
||||
@ -367,7 +367,6 @@ namespace Optimization {
|
||||
assert(0);
|
||||
}
|
||||
|
||||
|
||||
}
|
||||
|
||||
//////////////////////////////////////////////////////////////////////////////////////
|
||||
|
@ -78,11 +78,18 @@ namespace Grid {
|
||||
typedef typename RealPart < Scalar_type >::type Real;
|
||||
typedef Vector_type vector_type;
|
||||
typedef Scalar_type scalar_type;
|
||||
|
||||
|
||||
typedef union conv_t_union {
|
||||
Vector_type v;
|
||||
Scalar_type s[sizeof(Vector_type)/sizeof(Scalar_type)];
|
||||
conv_t_union(){};
|
||||
} conv_t;
|
||||
|
||||
|
||||
Vector_type v;
|
||||
|
||||
static inline int Nsimd(void) { return sizeof(Vector_type)/sizeof(Scalar_type);}
|
||||
|
||||
|
||||
Grid_simd& operator=(const Grid_simd&& rhs){v=rhs.v;return *this;};
|
||||
Grid_simd& operator=(const Grid_simd& rhs){v=rhs.v;return *this;}; //faster than not declaring it and leaving to the compiler
|
||||
@ -192,6 +199,27 @@ namespace Grid {
|
||||
return *this;
|
||||
}
|
||||
|
||||
|
||||
///////////////////////////////////////
|
||||
// Not all functions are supported
|
||||
// through SIMD and must breakout to
|
||||
// scalar type and back again. This
|
||||
// provides support
|
||||
///////////////////////////////////////
|
||||
|
||||
template<class functor> friend inline Grid_simd SimdApply (const functor &func,const Grid_simd &v) {
|
||||
|
||||
Grid_simd ret;
|
||||
Grid_simd::conv_t conv;
|
||||
|
||||
conv.v = v.v;
|
||||
for(int i=0;i<Nsimd();i++){
|
||||
conv.s[i]=func(conv.s[i]);
|
||||
}
|
||||
ret.v = conv.v;
|
||||
return ret;
|
||||
}
|
||||
|
||||
////////////////////////////////////////////////////////////////////
|
||||
// General permute; assumes vector length is same across
|
||||
// all subtypes; may not be a good assumption, but could
|
||||
|
60
lib/simd/Grid_vector_unops.h
Normal file
60
lib/simd/Grid_vector_unops.h
Normal file
@ -0,0 +1,60 @@
|
||||
#ifndef GRID_VECTOR_UNOPS
|
||||
#define GRID_VECTOR_UNOPS
|
||||
|
||||
namespace Grid {
|
||||
|
||||
template<class scalar> struct SqrtRealFunctor {
|
||||
scalar operator()(const scalar &a) const {
|
||||
return sqrt(real(a));
|
||||
}
|
||||
};
|
||||
|
||||
template<class scalar> struct RSqrtRealFunctor {
|
||||
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));
|
||||
}
|
||||
};
|
||||
|
||||
template<class scalar> struct SinRealFunctor {
|
||||
scalar operator()(const scalar &a) const {
|
||||
return sin(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);
|
||||
}
|
||||
};
|
||||
|
||||
template < class S, class V >
|
||||
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) {
|
||||
return SimdApply(RSqrtRealFunctor<S>(),r);
|
||||
}
|
||||
template < class S, class V >
|
||||
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) {
|
||||
return SimdApply(CosRealFunctor<S>(),r);
|
||||
}
|
||||
template < class S, class V >
|
||||
inline Grid_simd<S,V> pow(const Grid_simd<S,V> &r,double y) {
|
||||
return SimdApply(PowRealFunctor<S>(y),r);
|
||||
}
|
||||
|
||||
}
|
||||
#endif
|
Reference in New Issue
Block a user