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_ldlt(int rows, int cols)
00033 {
00034 static bool odd = true;
00035 odd = !odd;
00036 double density = std::max(8./(rows*cols), 0.01);
00037 typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
00038 typedef Matrix<Scalar,Dynamic,1> DenseVector;
00039
00040 SparseMatrix<Scalar> m2(rows, cols);
00041 DenseMatrix refMat2(rows, cols);
00042
00043 DenseVector b = DenseVector::Random(cols);
00044 DenseVector refX(cols), x(cols);
00045
00046 initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag|MakeUpperTriangular, 0, 0);
00047
00048 SparseMatrix<Scalar> m3 = m2 * m2.adjoint(), m3_lo(rows,rows), m3_up(rows,rows);
00049 DenseMatrix refMat3 = refMat2 * refMat2.adjoint();
00050
00051 refX = refMat3.template selfadjointView<Upper>().ldlt().solve(b);
00052 typedef SparseMatrix<Scalar,Upper|SelfAdjoint> SparseSelfAdjointMatrix;
00053 x = b;
00054 SparseLDLT<SparseSelfAdjointMatrix> ldlt(m3);
00055 if (ldlt.succeeded())
00056 ldlt.solveInPlace(x);
00057 else
00058 std::cerr << "warning LDLT failed\n";
00059
00060 VERIFY_IS_APPROX(refMat3.template selfadjointView<Upper>() * x, b);
00061 VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LDLT: default");
00062
00063 #ifdef EIGEN_CHOLMOD_SUPPORT
00064 {
00065 x = b;
00066 SparseLDLT<SparseSelfAdjointMatrix, Cholmod> ldlt2(m3);
00067 if (ldlt2.succeeded())
00068 {
00069 ldlt2.solveInPlace(x);
00070 VERIFY_IS_APPROX(refMat3.template selfadjointView<Upper>() * x, b);
00071 VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LDLT: cholmod solveInPlace");
00072
00073 x = ldlt2.solve(b);
00074 VERIFY_IS_APPROX(refMat3.template selfadjointView<Upper>() * x, b);
00075 VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LDLT: cholmod solve");
00076 }
00077 else
00078 std::cerr << "warning LDLT failed\n";
00079 }
00080 #endif
00081
00082
00083
00084
00085
00086 {
00087 SparseMatrix<Scalar> m2(rows, cols);
00088 DenseMatrix refMat2(rows, cols);
00089
00090 DenseVector b = DenseVector::Random(cols);
00091 DenseVector ref_x(cols), x(cols);
00092 DenseMatrix B = DenseMatrix::Random(rows,cols);
00093 DenseMatrix ref_X(rows,cols), X(rows,cols);
00094
00095 initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag|MakeLowerTriangular, 0, 0);
00096
00097 for(int i=0; i<rows; ++i)
00098 m2.coeffRef(i,i) = refMat2(i,i) = internal::abs(internal::real(refMat2(i,i)));
00099
00100
00101 SparseMatrix<Scalar> m3 = m2 * m2.adjoint(), m3_lo(rows,rows), m3_up(rows,rows);
00102 DenseMatrix refMat3 = refMat2 * refMat2.adjoint();
00103
00104 m3_lo.template selfadjointView<Lower>().rankUpdate(m2,0);
00105 m3_up.template selfadjointView<Upper>().rankUpdate(m2,0);
00106
00107
00108 ref_x = refMat3.template selfadjointView<Lower>().llt().solve(b);
00109
00110 x = SimplicialCholesky<SparseMatrix<Scalar>, Lower>().setMode(odd ? SimplicialCholeskyLLt : SimplicialCholeskyLDLt).compute(m3).solve(b);
00111 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "SimplicialCholesky: solve, full storage, lower, single dense rhs");
00112
00113 x = SimplicialCholesky<SparseMatrix<Scalar>, Upper>().setMode(odd ? SimplicialCholeskyLLt : SimplicialCholeskyLDLt).compute(m3).solve(b);
00114 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "SimplicialCholesky: solve, full storage, upper, single dense rhs");
00115
00116 x = SimplicialCholesky<SparseMatrix<Scalar>, Lower>(m3_lo).solve(b);
00117 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "SimplicialCholesky: solve, lower only, single dense rhs");
00118
00119 x = SimplicialCholesky<SparseMatrix<Scalar>, Upper>(m3_up).solve(b);
00120 VERIFY(ref_x.isApprox(x,test_precision<Scalar>()) && "SimplicialCholesky: solve, upper only, single dense rhs");
00121
00122
00123
00124 ref_X = refMat3.template selfadjointView<Lower>().llt().solve(B);
00125
00126 X = SimplicialCholesky<SparseMatrix<Scalar>, Lower>().setMode(odd ? SimplicialCholeskyLLt : SimplicialCholeskyLDLt).compute(m3).solve(B);
00127 VERIFY(ref_X.isApprox(X,test_precision<Scalar>()) && "SimplicialCholesky: solve, full storage, lower, multiple dense rhs");
00128
00129 X = SimplicialCholesky<SparseMatrix<Scalar>, Upper>().setMode(odd ? SimplicialCholeskyLLt : SimplicialCholeskyLDLt).compute(m3).solve(B);
00130 VERIFY(ref_X.isApprox(X,test_precision<Scalar>()) && "SimplicialCholesky: solve, full storage, upper, multiple dense rhs");
00131
00132
00133
00134
00135
00136
00137
00138
00139
00140
00141
00142
00143
00144
00145 }
00146
00147
00148
00149
00150
00151
00152
00153
00154
00155
00156
00157
00158
00159
00160
00161
00162
00163
00164
00165 }
00166
00167 void test_sparse_ldlt()
00168 {
00169 for(int i = 0; i < g_repeat; i++) {
00170 CALL_SUBTEST_1(sparse_ldlt<double>(8, 8) );
00171 int s = internal::random<int>(1,300);
00172 CALL_SUBTEST_2(sparse_ldlt<std::complex<double> >(s,s) );
00173 CALL_SUBTEST_1(sparse_ldlt<double>(s,s) );
00174 }
00175 }