aboutsummaryrefslogtreecommitdiffhomepage
path: root/unsupported/Eigen/MPRealSupport
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2014-07-17 12:00:56 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2014-07-17 12:00:56 +0200
commit84ad8ce7e30e26eb235906cd569d1f2b7f3798c5 (patch)
treef4463145f6fe589ea157ca775d3f1793cd2e13aa /unsupported/Eigen/MPRealSupport
parent40b74411e4ba49990bef69175a866f53587f8c59 (diff)
Fix bug #770: workaround thread safety in mpreal
Diffstat (limited to 'unsupported/Eigen/MPRealSupport')
-rw-r--r--unsupported/Eigen/MPRealSupport28
1 files changed, 14 insertions, 14 deletions
diff --git a/unsupported/Eigen/MPRealSupport b/unsupported/Eigen/MPRealSupport
index 6a65a601b..de4f11934 100644
--- a/unsupported/Eigen/MPRealSupport
+++ b/unsupported/Eigen/MPRealSupport
@@ -146,8 +146,8 @@ int main()
};
};
- template<typename Index, int mr, int nr, bool ConjugateLhs, bool ConjugateRhs>
- struct gebp_kernel<mpfr::mpreal,mpfr::mpreal,Index,mr,nr,ConjugateLhs,ConjugateRhs>
+ template<typename Index, bool ConjugateLhs, bool ConjugateRhs>
+ struct gebp_kernel<mpfr::mpreal,mpfr::mpreal,Index,1,1,ConjugateLhs,ConjugateRhs>
{
typedef mpfr::mpreal mpreal;
@@ -155,34 +155,34 @@ int main()
void operator()(mpreal* res, Index resStride, const mpreal* blockA, const mpreal* blockB, Index rows, Index depth, Index cols, mpreal alpha,
Index strideA=-1, Index strideB=-1, Index offsetA=0, Index offsetB=0)
{
- mpreal acc1, tmp;
+ mpreal acc1(0,mpfr_get_prec(blockA[0].mpfr_srcptr())),
+ tmp (0,mpfr_get_prec(blockA[0].mpfr_srcptr()));
if(strideA==-1) strideA = depth;
if(strideB==-1) strideB = depth;
- for(Index j=0; j<cols; j+=nr)
+ for(Index i=0; i<rows; ++i)
{
- Index actual_nr = (std::min<Index>)(nr,cols-j);
- mpreal *C1 = res + j*resStride;
- for(Index i=0; i<rows; i++)
+ for(Index j=0; j<cols; ++j)
{
- mpreal *B = const_cast<mpreal*>(blockB) + j*strideB + offsetB*actual_nr;
- mpreal *A = const_cast<mpreal*>(blockA) + i*strideA + offsetA;
+ mpreal *C1 = res + j*resStride;
+
+ const mpreal *A = blockA + i*strideA + offsetA;
+ const mpreal *B = blockB + j*strideB + offsetB;
+
acc1 = 0;
for(Index k=0; k<depth; k++)
{
- mpfr_mul(tmp.mpfr_ptr(), A[k].mpfr_ptr(), B[0].mpfr_ptr(), mpreal::get_default_rnd());
+ mpfr_mul(tmp.mpfr_ptr(), A[k].mpfr_srcptr(), B[k].mpfr_srcptr(), mpreal::get_default_rnd());
mpfr_add(acc1.mpfr_ptr(), acc1.mpfr_ptr(), tmp.mpfr_ptr(), mpreal::get_default_rnd());
- B+=actual_nr;
}
- mpfr_mul(acc1.mpfr_ptr(), acc1.mpfr_ptr(), alpha.mpfr_ptr(), mpreal::get_default_rnd());
- mpfr_add(C1[i].mpfr_ptr(), C1[i].mpfr_ptr(), acc1.mpfr_ptr(), mpreal::get_default_rnd());
+ 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());
}
}
}
};
-
} // end namespace internal
}