mirror of
https://gitlab.com/libeigen/eigen.git
synced 2025-04-24 02:29:33 +08:00
Extend mpreal unit test to check LLT with complexes.
This commit is contained in:
parent
a354c3ca59
commit
55b4fd1d40
@ -12,6 +12,7 @@ void test_mpreal_support()
|
||||
// set precision to 256 bits (double has only 53 bits)
|
||||
mpreal::set_default_prec(256);
|
||||
typedef Matrix<mpreal,Eigen::Dynamic,Eigen::Dynamic> MatrixXmp;
|
||||
typedef Matrix<std::complex<mpreal>,Eigen::Dynamic,Eigen::Dynamic> MatrixXcmp;
|
||||
|
||||
std::cerr << "epsilon = " << NumTraits<mpreal>::epsilon() << "\n";
|
||||
std::cerr << "dummy_precision = " << NumTraits<mpreal>::dummy_precision() << "\n";
|
||||
@ -25,6 +26,10 @@ void test_mpreal_support()
|
||||
MatrixXmp B = MatrixXmp::Random(s,s);
|
||||
MatrixXmp S = A.adjoint() * A;
|
||||
MatrixXmp X;
|
||||
MatrixXcmp Ac = MatrixXcmp::Random(s,s);
|
||||
MatrixXcmp Bc = MatrixXcmp::Random(s,s);
|
||||
MatrixXcmp Sc = Ac.adjoint() * Ac;
|
||||
MatrixXcmp Xc;
|
||||
|
||||
// Basic stuffs
|
||||
VERIFY_IS_APPROX(A.real(), A);
|
||||
@ -37,6 +42,9 @@ void test_mpreal_support()
|
||||
// Cholesky
|
||||
X = S.selfadjointView<Lower>().llt().solve(B);
|
||||
VERIFY_IS_APPROX((S.selfadjointView<Lower>()*X).eval(),B);
|
||||
|
||||
Xc = Sc.selfadjointView<Lower>().llt().solve(Bc);
|
||||
VERIFY_IS_APPROX((Sc.selfadjointView<Lower>()*Xc).eval(),Bc);
|
||||
|
||||
// partial LU
|
||||
X = A.lu().solve(B);
|
||||
|
Loading…
x
Reference in New Issue
Block a user