aboutsummaryrefslogtreecommitdiffhomepage
diff options
context:
space:
mode:
authorGravatar Gael Guennebaud <g.gael@free.fr>2012-06-21 10:02:32 +0200
committerGravatar Gael Guennebaud <g.gael@free.fr>2012-06-21 10:02:32 +0200
commit5b5f3ecafa25aee1af27ea3503ad9e0e1ee49d04 (patch)
tree055a0f8199fb5a526b46d8bb4be2a4c1c5a1d7b9
parent7380592bc23aa9d331b99bac22bbb5fb4c854639 (diff)
MPreal: extended unit test, remove useless internal overloads, add support for internal::cast (needed for printing)
-rw-r--r--unsupported/Eigen/MPRealSupport27
-rw-r--r--unsupported/test/mpreal_support.cpp17
2 files changed, 29 insertions, 15 deletions
diff --git a/unsupported/Eigen/MPRealSupport b/unsupported/Eigen/MPRealSupport
index 44d75856e..b5fe5e404 100644
--- a/unsupported/Eigen/MPRealSupport
+++ b/unsupported/Eigen/MPRealSupport
@@ -107,7 +107,7 @@ int main()
}
};
- namespace internal {
+namespace internal {
template<> mpfr::mpreal random<mpfr::mpreal>()
{
@@ -133,18 +133,6 @@ int main()
return a + (b-a) * random<mpfr::mpreal>();
}
- template<> struct conj_impl<mpfr::mpreal> { inline static const mpfr::mpreal& run(const mpfr::mpreal& x) { return x; } };
- template<> struct real_impl<mpfr::mpreal> { inline static const mpfr::mpreal& run(const mpfr::mpreal& x) { return x; } };
- template<> struct imag_impl<mpfr::mpreal> { inline static const mpfr::mpreal run(const mpfr::mpreal&) { return mpfr::mpreal(0); } };
- template<> struct abs_impl<mpfr::mpreal> { inline static const mpfr::mpreal run(const mpfr::mpreal& x) { return mpfr::fabs(x); } };
- template<> struct abs2_impl<mpfr::mpreal> { inline static const mpfr::mpreal run(const mpfr::mpreal& x) { return x*x; } };
- template<> struct sqrt_impl<mpfr::mpreal> { inline static const mpfr::mpreal run(const mpfr::mpreal& x) { return mpfr::sqrt(x); } };
- template<> struct exp_impl<mpfr::mpreal> { inline static const mpfr::mpreal run(const mpfr::mpreal& x) { return mpfr::exp(x); } };
- template<> struct log_impl<mpfr::mpreal> { inline static const mpfr::mpreal run(const mpfr::mpreal& x) { return mpfr::log(x); } };
- template<> struct sin_impl<mpfr::mpreal> { inline static const mpfr::mpreal run(const mpfr::mpreal& x) { return mpfr::sin(x); } };
- template<> struct cos_impl<mpfr::mpreal> { inline static const mpfr::mpreal run(const mpfr::mpreal& x) { return mpfr::cos(x); } };
- template<> struct pow_impl<mpfr::mpreal> { inline static const mpfr::mpreal run(const mpfr::mpreal& x, const mpfr::mpreal& y) { return mpfr::pow(x, y); } };
-
bool isMuchSmallerThan(const mpfr::mpreal& a, const mpfr::mpreal& b, const mpfr::mpreal& prec)
{
return mpfr::abs(a) <= mpfr::abs(b) * prec;
@@ -159,8 +147,17 @@ int main()
{
return a <= b || isApprox(a, b, prec);
}
-
- } // end namespace internal
+
+ template<> inline long double cast<mpfr::mpreal,long double>(const mpfr::mpreal& x)
+ { return x.toLDouble(); }
+ template<> inline double cast<mpfr::mpreal,double>(const mpfr::mpreal& x)
+ { return x.toDouble(); }
+ template<> inline long cast<mpfr::mpreal,long>(const mpfr::mpreal& x)
+ { return x.toLong(); }
+ template<> inline int cast<mpfr::mpreal,int>(const mpfr::mpreal& x)
+ { return int(x.toLong()); }
+
+} // end namespace internal
}
#endif // EIGEN_MPREALSUPPORT_MODULE_H
diff --git a/unsupported/test/mpreal_support.cpp b/unsupported/test/mpreal_support.cpp
index 53d388821..551af9db8 100644
--- a/unsupported/test/mpreal_support.cpp
+++ b/unsupported/test/mpreal_support.cpp
@@ -2,6 +2,7 @@
#include <Eigen/MPRealSupport>
#include <Eigen/LU>
#include <Eigen/Eigenvalues>
+#include <sstream>
using namespace mpfr;
using namespace std;
@@ -24,6 +25,15 @@ void test_mpreal_support()
MatrixXmp B = MatrixXmp::Random(s,s);
MatrixXmp S = A.adjoint() * A;
MatrixXmp X;
+
+ // Basic stuffs
+ VERIFY_IS_APPROX(A.real(), A);
+ VERIFY(Eigen::internal::isApprox(A.array().abs2().sum(), A.squaredNorm()));
+ VERIFY_IS_APPROX(A.array().exp(), exp(A.array()));
+ VERIFY_IS_APPROX(A.array().abs2().sqrt(), A.array().abs());
+ VERIFY_IS_APPROX(A.array().sin(), sin(A.array()));
+ VERIFY_IS_APPROX(A.array().cos(), cos(A.array()));
+
// Cholesky
X = S.selfadjointView<Lower>().llt().solve(B);
@@ -39,6 +49,13 @@ void test_mpreal_support()
VERIFY_IS_APPROX((S.selfadjointView<Lower>() * eig.eigenvectors()),
eig.eigenvectors() * eig.eigenvalues().asDiagonal());
}
+
+ {
+ MatrixXmp A(8,3); A.setRandom();
+ // test output (interesting things happen in this code)
+ std::stringstream stream;
+ stream << A;
+ }
}
extern "C" {