00001
00002
00003
00004
00005
00006
00007
00008
00009
00010
00011
00012
00013
00014
00015
00016
00017
00018
00019
00020
00021
00022
00023
00024
00025 #include "sparse.h"
00026 #include <Eigen/SparseExtra>
00027
00028 #ifdef EIGEN_CHOLMOD_SUPPORT
00029 #include <Eigen/CholmodSupport>
00030 #endif
00031
00032 template<typename Scalar> void sparse_llt(int rows, int cols)
00033 {
00034 double density = std::max(8./(rows*cols), 0.01);
00035 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
00036 typedef Matrix<Scalar,Dynamic,1> DenseVector;
00037
00038
00039 SparseMatrix<Scalar> m2(rows, cols);
00040 DenseMatrix refMat2(rows, cols);
00041
00042 DenseVector b = DenseVector::Random(cols);
00043 DenseVector ref_x(cols), x(cols);
00044 DenseMatrix B = DenseMatrix::Random(rows,cols);
00045 DenseMatrix ref_X(rows,cols), X(rows,cols);
00046
00047 initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag|MakeLowerTriangular, 0, 0);
00048
00049 for(int i=0; i<rows; ++i)
00050 m2.coeffRef(i,i) = refMat2(i,i) = internal::abs(internal::real(refMat2(i,i)));
00051
00052 ref_x = refMat2.template selfadjointView<Lower>().llt().solve(b);
00053 if (!NumTraits<Scalar>::IsComplex)
00054 {
00055 x = b;
00056 SparseLLT<SparseMatrix<Scalar> > (m2).solveInPlace(x);
00057 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "LLT: default");
00058 }
00059
00060 #ifdef EIGEN_CHOLMOD_SUPPORT
00061
00062 {
00063
00064 SparseMatrix<Scalar> m3 = m2.adjoint()*m2;
00065 DenseMatrix refMat3 = refMat2.adjoint()*refMat2;
00066
00067 ref_x = refMat3.template selfadjointView<Lower>().llt().solve(b);
00068
00069 x = b;
00070 SparseLLT<SparseMatrix<Scalar>, Cholmod>(m3).solveInPlace(x);
00071 VERIFY((m3*x).isApprox(b,test_precision<Scalar>()) && "LLT legacy: cholmod solveInPlace");
00072
00073 x = SparseLLT<SparseMatrix<Scalar>, Cholmod>(m3).solve(b);
00074 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "LLT legacy: cholmod solve");
00075 }
00076
00077
00078 {
00079
00080 SparseMatrix<Scalar> m3 = m2 * m2.adjoint(), m3_lo(rows,rows), m3_up(rows,rows);
00081 DenseMatrix refMat3 = refMat2 * refMat2.adjoint();
00082
00083 m3_lo.template selfadjointView<Lower>().rankUpdate(m2,0);
00084 m3_up.template selfadjointView<Upper>().rankUpdate(m2,0);
00085
00086
00087 ref_x = refMat3.template selfadjointView<Lower>().llt().solve(b);
00088
00089 x = CholmodDecomposition<SparseMatrix<Scalar>, Lower>(m3).solve(b);
00090 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "LLT: cholmod solve, single dense rhs");
00091
00092 x = CholmodDecomposition<SparseMatrix<Scalar>, Upper>(m3).solve(b);
00093 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "LLT: cholmod solve, single dense rhs");
00094
00095 x = CholmodDecomposition<SparseMatrix<Scalar>, Lower>(m3_lo).solve(b);
00096 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "LLT: cholmod solve, single dense rhs");
00097
00098 x = CholmodDecomposition<SparseMatrix<Scalar>, Upper>(m3_up).solve(b);
00099 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "LLT: cholmod solve, single dense rhs");
00100
00101
00102
00103 ref_X = refMat3.template selfadjointView<Lower>().llt().solve(B);
00104
00105 #ifndef EIGEN_DEFAULT_TO_ROW_MAJOR
00106
00107 X = CholmodDecomposition<SparseMatrix<Scalar>, Lower>(m3).solve(B);
00108 VERIFY(ref_X.isApprox(X,test_precision<Scalar>()) && "LLT: cholmod solve, multiple dense rhs");
00109
00110 X = CholmodDecomposition<SparseMatrix<Scalar>, Upper>(m3).solve(B);
00111 VERIFY(ref_X.isApprox(X,test_precision<Scalar>()) && "LLT: cholmod solve, multiple dense rhs");
00112 #endif
00113
00114
00115
00116 SparseMatrix<Scalar> spB(rows,cols), spX(rows,cols);
00117 B.diagonal().array() += 1;
00118 spB = B.sparseView(0.5,1);
00119
00120 ref_X = refMat3.template selfadjointView<Lower>().llt().solve(DenseMatrix(spB));
00121
00122 spX = CholmodDecomposition<SparseMatrix<Scalar>, Lower>(m3).solve(spB);
00123 VERIFY(ref_X.isApprox(spX.toDense(),test_precision<Scalar>()) && "LLT: cholmod solve, multiple sparse rhs");
00124
00125 spX = CholmodDecomposition<SparseMatrix<Scalar>, Upper>(m3).solve(spB);
00126 VERIFY(ref_X.isApprox(spX.toDense(),test_precision<Scalar>()) && "LLT: cholmod solve, multiple sparse rhs");
00127 }
00128 #endif
00129
00130 }
00131
00132 void test_sparse_llt()
00133 {
00134 for(int i = 0; i < g_repeat; i++) {
00135 CALL_SUBTEST_1(sparse_llt<double>(8, 8) );
00136 int s = internal::random<int>(1,300);
00137 CALL_SUBTEST_2(sparse_llt<std::complex<double> >(s,s) );
00138 CALL_SUBTEST_1(sparse_llt<double>(s,s) );
00139 }
00140 }