aboutsummaryrefslogtreecommitdiffhomepage
path: root/unsupported/Eigen/MPRealSupport
diff options
context:
space:
mode:
authorGravatar Christoph Hertzberg <chtz@informatik.uni-bremen.de>2015-02-28 16:41:00 +0100
committerGravatar Christoph Hertzberg <chtz@informatik.uni-bremen.de>2015-02-28 16:41:00 +0100
commit682196e9fcdef7fb329fde833e30e91bb3f89077 (patch)
treeb478f6363c6283767aaec87f0943bb2bb950e62f /unsupported/Eigen/MPRealSupport
parent33f40b2883630b5ddec3bab5730691b40120faa6 (diff)
Fixed MPRealSupport
Diffstat (limited to 'unsupported/Eigen/MPRealSupport')
-rw-r--r--unsupported/Eigen/MPRealSupport22
1 files changed, 16 insertions, 6 deletions
diff --git a/unsupported/Eigen/MPRealSupport b/unsupported/Eigen/MPRealSupport
index 8e42965a3..89036886b 100644
--- a/unsupported/Eigen/MPRealSupport
+++ b/unsupported/Eigen/MPRealSupport
@@ -141,20 +141,32 @@ int main()
public:
typedef mpfr::mpreal ResScalar;
enum {
+ Vectorizable = false,
+ LhsPacketSize = 1,
+ RhsPacketSize = 1,
+ ResPacketSize = 1,
+ NumberOfRegisters = 1,
nr = 1,
mr = 1,
LhsProgress = 1,
RhsProgress = 1
};
+ typedef ResScalar LhsPacket;
+ typedef ResScalar RhsPacket;
+ typedef ResScalar ResPacket;
+
};
- template<typename Index, bool ConjugateLhs, bool ConjugateRhs>
- struct gebp_kernel<mpfr::mpreal,mpfr::mpreal,Index,1,1,ConjugateLhs,ConjugateRhs>
+
+
+ template<typename Index, typename DataMapper, bool ConjugateLhs, bool ConjugateRhs>
+ struct gebp_kernel<mpfr::mpreal,mpfr::mpreal,Index,DataMapper,1,1,ConjugateLhs,ConjugateRhs>
{
typedef mpfr::mpreal mpreal;
EIGEN_DONT_INLINE
- void operator()(mpreal* res, Index resStride, const mpreal* blockA, const mpreal* blockB, Index rows, Index depth, Index cols, mpreal alpha,
+ void operator()(const DataMapper& res, const mpreal* blockA, const mpreal* blockB,
+ Index rows, Index depth, Index cols, const mpreal& alpha,
Index strideA=-1, Index strideB=-1, Index offsetA=0, Index offsetB=0)
{
if(rows==0 || cols==0 || depth==0)
@@ -170,8 +182,6 @@ int main()
{
for(Index j=0; j<cols; ++j)
{
- mpreal *C1 = res + j*resStride;
-
const mpreal *A = blockA + i*strideA + offsetA;
const mpreal *B = blockB + j*strideB + offsetB;
@@ -183,7 +193,7 @@ int main()
}
mpfr_mul(acc1.mpfr_ptr(), acc1.mpfr_srcptr(), alpha.mpfr_srcptr(), mpreal::get_default_rnd());
- mpfr_add(C1[i].mpfr_ptr(), C1[i].mpfr_srcptr(), acc1.mpfr_srcptr(), mpreal::get_default_rnd());
+ mpfr_add(res(i,j).mpfr_ptr(), res(i,j).mpfr_srcptr(), acc1.mpfr_srcptr(), mpreal::get_default_rnd());
}
}
}