xref: /aosp_15_r20/external/eigen/test/cholesky.cpp (revision bf2c37156dfe67e5dfebd6d394bad8b2ab5804d4)
1*bf2c3715SXin Li // This file is part of Eigen, a lightweight C++ template library
2*bf2c3715SXin Li // for linear algebra.
3*bf2c3715SXin Li //
4*bf2c3715SXin Li // Copyright (C) 2008 Gael Guennebaud <[email protected]>
5*bf2c3715SXin Li //
6*bf2c3715SXin Li // This Source Code Form is subject to the terms of the Mozilla
7*bf2c3715SXin Li // Public License v. 2.0. If a copy of the MPL was not distributed
8*bf2c3715SXin Li // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9*bf2c3715SXin Li 
10*bf2c3715SXin Li #define TEST_ENABLE_TEMPORARY_TRACKING
11*bf2c3715SXin Li 
12*bf2c3715SXin Li #include "main.h"
13*bf2c3715SXin Li #include <Eigen/Cholesky>
14*bf2c3715SXin Li #include <Eigen/QR>
15*bf2c3715SXin Li #include "solverbase.h"
16*bf2c3715SXin Li 
17*bf2c3715SXin Li template<typename MatrixType, int UpLo>
matrix_l1_norm(const MatrixType & m)18*bf2c3715SXin Li typename MatrixType::RealScalar matrix_l1_norm(const MatrixType& m) {
19*bf2c3715SXin Li   if(m.cols()==0) return typename MatrixType::RealScalar(0);
20*bf2c3715SXin Li   MatrixType symm = m.template selfadjointView<UpLo>();
21*bf2c3715SXin Li   return symm.cwiseAbs().colwise().sum().maxCoeff();
22*bf2c3715SXin Li }
23*bf2c3715SXin Li 
test_chol_update(const MatrixType & symm)24*bf2c3715SXin Li template<typename MatrixType,template <typename,int> class CholType> void test_chol_update(const MatrixType& symm)
25*bf2c3715SXin Li {
26*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
27*bf2c3715SXin Li   typedef typename MatrixType::RealScalar RealScalar;
28*bf2c3715SXin Li   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
29*bf2c3715SXin Li 
30*bf2c3715SXin Li   MatrixType symmLo = symm.template triangularView<Lower>();
31*bf2c3715SXin Li   MatrixType symmUp = symm.template triangularView<Upper>();
32*bf2c3715SXin Li   MatrixType symmCpy = symm;
33*bf2c3715SXin Li 
34*bf2c3715SXin Li   CholType<MatrixType,Lower> chollo(symmLo);
35*bf2c3715SXin Li   CholType<MatrixType,Upper> cholup(symmUp);
36*bf2c3715SXin Li 
37*bf2c3715SXin Li   for (int k=0; k<10; ++k)
38*bf2c3715SXin Li   {
39*bf2c3715SXin Li     VectorType vec = VectorType::Random(symm.rows());
40*bf2c3715SXin Li     RealScalar sigma = internal::random<RealScalar>();
41*bf2c3715SXin Li     symmCpy += sigma * vec * vec.adjoint();
42*bf2c3715SXin Li 
43*bf2c3715SXin Li     // we are doing some downdates, so it might be the case that the matrix is not SPD anymore
44*bf2c3715SXin Li     CholType<MatrixType,Lower> chol(symmCpy);
45*bf2c3715SXin Li     if(chol.info()!=Success)
46*bf2c3715SXin Li       break;
47*bf2c3715SXin Li 
48*bf2c3715SXin Li     chollo.rankUpdate(vec, sigma);
49*bf2c3715SXin Li     VERIFY_IS_APPROX(symmCpy, chollo.reconstructedMatrix());
50*bf2c3715SXin Li 
51*bf2c3715SXin Li     cholup.rankUpdate(vec, sigma);
52*bf2c3715SXin Li     VERIFY_IS_APPROX(symmCpy, cholup.reconstructedMatrix());
53*bf2c3715SXin Li   }
54*bf2c3715SXin Li }
55*bf2c3715SXin Li 
cholesky(const MatrixType & m)56*bf2c3715SXin Li template<typename MatrixType> void cholesky(const MatrixType& m)
57*bf2c3715SXin Li {
58*bf2c3715SXin Li   /* this test covers the following files:
59*bf2c3715SXin Li      LLT.h LDLT.h
60*bf2c3715SXin Li   */
61*bf2c3715SXin Li   Index rows = m.rows();
62*bf2c3715SXin Li   Index cols = m.cols();
63*bf2c3715SXin Li 
64*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
65*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
66*bf2c3715SXin Li   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> SquareMatrixType;
67*bf2c3715SXin Li   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
68*bf2c3715SXin Li 
69*bf2c3715SXin Li   MatrixType a0 = MatrixType::Random(rows,cols);
70*bf2c3715SXin Li   VectorType vecB = VectorType::Random(rows), vecX(rows);
71*bf2c3715SXin Li   MatrixType matB = MatrixType::Random(rows,cols), matX(rows,cols);
72*bf2c3715SXin Li   SquareMatrixType symm =  a0 * a0.adjoint();
73*bf2c3715SXin Li   // let's make sure the matrix is not singular or near singular
74*bf2c3715SXin Li   for (int k=0; k<3; ++k)
75*bf2c3715SXin Li   {
76*bf2c3715SXin Li     MatrixType a1 = MatrixType::Random(rows,cols);
77*bf2c3715SXin Li     symm += a1 * a1.adjoint();
78*bf2c3715SXin Li   }
79*bf2c3715SXin Li 
80*bf2c3715SXin Li   {
81*bf2c3715SXin Li     STATIC_CHECK(( internal::is_same<typename LLT<MatrixType,Lower>::StorageIndex,int>::value ));
82*bf2c3715SXin Li     STATIC_CHECK(( internal::is_same<typename LLT<MatrixType,Upper>::StorageIndex,int>::value ));
83*bf2c3715SXin Li 
84*bf2c3715SXin Li     SquareMatrixType symmUp = symm.template triangularView<Upper>();
85*bf2c3715SXin Li     SquareMatrixType symmLo = symm.template triangularView<Lower>();
86*bf2c3715SXin Li 
87*bf2c3715SXin Li     LLT<SquareMatrixType,Lower> chollo(symmLo);
88*bf2c3715SXin Li     VERIFY_IS_APPROX(symm, chollo.reconstructedMatrix());
89*bf2c3715SXin Li 
90*bf2c3715SXin Li     check_solverbase<VectorType, VectorType>(symm, chollo, rows, rows, 1);
91*bf2c3715SXin Li     check_solverbase<MatrixType, MatrixType>(symm, chollo, rows, cols, rows);
92*bf2c3715SXin Li 
93*bf2c3715SXin Li     const MatrixType symmLo_inverse = chollo.solve(MatrixType::Identity(rows,cols));
94*bf2c3715SXin Li     RealScalar rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Lower>(symmLo)) /
95*bf2c3715SXin Li                              matrix_l1_norm<MatrixType, Lower>(symmLo_inverse);
96*bf2c3715SXin Li     RealScalar rcond_est = chollo.rcond();
97*bf2c3715SXin Li     // Verify that the estimated condition number is within a factor of 10 of the
98*bf2c3715SXin Li     // truth.
99*bf2c3715SXin Li     VERIFY(rcond_est >= rcond / 10 && rcond_est <= rcond * 10);
100*bf2c3715SXin Li 
101*bf2c3715SXin Li     // test the upper mode
102*bf2c3715SXin Li     LLT<SquareMatrixType,Upper> cholup(symmUp);
103*bf2c3715SXin Li     VERIFY_IS_APPROX(symm, cholup.reconstructedMatrix());
104*bf2c3715SXin Li     vecX = cholup.solve(vecB);
105*bf2c3715SXin Li     VERIFY_IS_APPROX(symm * vecX, vecB);
106*bf2c3715SXin Li     matX = cholup.solve(matB);
107*bf2c3715SXin Li     VERIFY_IS_APPROX(symm * matX, matB);
108*bf2c3715SXin Li 
109*bf2c3715SXin Li     // Verify that the estimated condition number is within a factor of 10 of the
110*bf2c3715SXin Li     // truth.
111*bf2c3715SXin Li     const MatrixType symmUp_inverse = cholup.solve(MatrixType::Identity(rows,cols));
112*bf2c3715SXin Li     rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Upper>(symmUp)) /
113*bf2c3715SXin Li                              matrix_l1_norm<MatrixType, Upper>(symmUp_inverse);
114*bf2c3715SXin Li     rcond_est = cholup.rcond();
115*bf2c3715SXin Li     VERIFY(rcond_est >= rcond / 10 && rcond_est <= rcond * 10);
116*bf2c3715SXin Li 
117*bf2c3715SXin Li 
118*bf2c3715SXin Li     MatrixType neg = -symmLo;
119*bf2c3715SXin Li     chollo.compute(neg);
120*bf2c3715SXin Li     VERIFY(neg.size()==0 || chollo.info()==NumericalIssue);
121*bf2c3715SXin Li 
122*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(chollo.matrixL().transpose().conjugate()), MatrixType(chollo.matrixU()));
123*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(chollo.matrixU().transpose().conjugate()), MatrixType(chollo.matrixL()));
124*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(cholup.matrixL().transpose().conjugate()), MatrixType(cholup.matrixU()));
125*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(cholup.matrixU().transpose().conjugate()), MatrixType(cholup.matrixL()));
126*bf2c3715SXin Li 
127*bf2c3715SXin Li     // test some special use cases of SelfCwiseBinaryOp:
128*bf2c3715SXin Li     MatrixType m1 = MatrixType::Random(rows,cols), m2(rows,cols);
129*bf2c3715SXin Li     m2 = m1;
130*bf2c3715SXin Li     m2 += symmLo.template selfadjointView<Lower>().llt().solve(matB);
131*bf2c3715SXin Li     VERIFY_IS_APPROX(m2, m1 + symmLo.template selfadjointView<Lower>().llt().solve(matB));
132*bf2c3715SXin Li     m2 = m1;
133*bf2c3715SXin Li     m2 -= symmLo.template selfadjointView<Lower>().llt().solve(matB);
134*bf2c3715SXin Li     VERIFY_IS_APPROX(m2, m1 - symmLo.template selfadjointView<Lower>().llt().solve(matB));
135*bf2c3715SXin Li     m2 = m1;
136*bf2c3715SXin Li     m2.noalias() += symmLo.template selfadjointView<Lower>().llt().solve(matB);
137*bf2c3715SXin Li     VERIFY_IS_APPROX(m2, m1 + symmLo.template selfadjointView<Lower>().llt().solve(matB));
138*bf2c3715SXin Li     m2 = m1;
139*bf2c3715SXin Li     m2.noalias() -= symmLo.template selfadjointView<Lower>().llt().solve(matB);
140*bf2c3715SXin Li     VERIFY_IS_APPROX(m2, m1 - symmLo.template selfadjointView<Lower>().llt().solve(matB));
141*bf2c3715SXin Li   }
142*bf2c3715SXin Li 
143*bf2c3715SXin Li   // LDLT
144*bf2c3715SXin Li   {
145*bf2c3715SXin Li     STATIC_CHECK(( internal::is_same<typename LDLT<MatrixType,Lower>::StorageIndex,int>::value ));
146*bf2c3715SXin Li     STATIC_CHECK(( internal::is_same<typename LDLT<MatrixType,Upper>::StorageIndex,int>::value ));
147*bf2c3715SXin Li 
148*bf2c3715SXin Li     int sign = internal::random<int>()%2 ? 1 : -1;
149*bf2c3715SXin Li 
150*bf2c3715SXin Li     if(sign == -1)
151*bf2c3715SXin Li     {
152*bf2c3715SXin Li       symm = -symm; // test a negative matrix
153*bf2c3715SXin Li     }
154*bf2c3715SXin Li 
155*bf2c3715SXin Li     SquareMatrixType symmUp = symm.template triangularView<Upper>();
156*bf2c3715SXin Li     SquareMatrixType symmLo = symm.template triangularView<Lower>();
157*bf2c3715SXin Li 
158*bf2c3715SXin Li     LDLT<SquareMatrixType,Lower> ldltlo(symmLo);
159*bf2c3715SXin Li     VERIFY(ldltlo.info()==Success);
160*bf2c3715SXin Li     VERIFY_IS_APPROX(symm, ldltlo.reconstructedMatrix());
161*bf2c3715SXin Li 
162*bf2c3715SXin Li     check_solverbase<VectorType, VectorType>(symm, ldltlo, rows, rows, 1);
163*bf2c3715SXin Li     check_solverbase<MatrixType, MatrixType>(symm, ldltlo, rows, cols, rows);
164*bf2c3715SXin Li 
165*bf2c3715SXin Li     const MatrixType symmLo_inverse = ldltlo.solve(MatrixType::Identity(rows,cols));
166*bf2c3715SXin Li     RealScalar rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Lower>(symmLo)) /
167*bf2c3715SXin Li                              matrix_l1_norm<MatrixType, Lower>(symmLo_inverse);
168*bf2c3715SXin Li     RealScalar rcond_est = ldltlo.rcond();
169*bf2c3715SXin Li     // Verify that the estimated condition number is within a factor of 10 of the
170*bf2c3715SXin Li     // truth.
171*bf2c3715SXin Li     VERIFY(rcond_est >= rcond / 10 && rcond_est <= rcond * 10);
172*bf2c3715SXin Li 
173*bf2c3715SXin Li 
174*bf2c3715SXin Li     LDLT<SquareMatrixType,Upper> ldltup(symmUp);
175*bf2c3715SXin Li     VERIFY(ldltup.info()==Success);
176*bf2c3715SXin Li     VERIFY_IS_APPROX(symm, ldltup.reconstructedMatrix());
177*bf2c3715SXin Li     vecX = ldltup.solve(vecB);
178*bf2c3715SXin Li     VERIFY_IS_APPROX(symm * vecX, vecB);
179*bf2c3715SXin Li     matX = ldltup.solve(matB);
180*bf2c3715SXin Li     VERIFY_IS_APPROX(symm * matX, matB);
181*bf2c3715SXin Li 
182*bf2c3715SXin Li     // Verify that the estimated condition number is within a factor of 10 of the
183*bf2c3715SXin Li     // truth.
184*bf2c3715SXin Li     const MatrixType symmUp_inverse = ldltup.solve(MatrixType::Identity(rows,cols));
185*bf2c3715SXin Li     rcond = (RealScalar(1) / matrix_l1_norm<MatrixType, Upper>(symmUp)) /
186*bf2c3715SXin Li                              matrix_l1_norm<MatrixType, Upper>(symmUp_inverse);
187*bf2c3715SXin Li     rcond_est = ldltup.rcond();
188*bf2c3715SXin Li     VERIFY(rcond_est >= rcond / 10 && rcond_est <= rcond * 10);
189*bf2c3715SXin Li 
190*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(ldltlo.matrixL().transpose().conjugate()), MatrixType(ldltlo.matrixU()));
191*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(ldltlo.matrixU().transpose().conjugate()), MatrixType(ldltlo.matrixL()));
192*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(ldltup.matrixL().transpose().conjugate()), MatrixType(ldltup.matrixU()));
193*bf2c3715SXin Li     VERIFY_IS_APPROX(MatrixType(ldltup.matrixU().transpose().conjugate()), MatrixType(ldltup.matrixL()));
194*bf2c3715SXin Li 
195*bf2c3715SXin Li     if(MatrixType::RowsAtCompileTime==Dynamic)
196*bf2c3715SXin Li     {
197*bf2c3715SXin Li       // note : each inplace permutation requires a small temporary vector (mask)
198*bf2c3715SXin Li 
199*bf2c3715SXin Li       // check inplace solve
200*bf2c3715SXin Li       matX = matB;
201*bf2c3715SXin Li       VERIFY_EVALUATION_COUNT(matX = ldltlo.solve(matX), 0);
202*bf2c3715SXin Li       VERIFY_IS_APPROX(matX, ldltlo.solve(matB).eval());
203*bf2c3715SXin Li 
204*bf2c3715SXin Li 
205*bf2c3715SXin Li       matX = matB;
206*bf2c3715SXin Li       VERIFY_EVALUATION_COUNT(matX = ldltup.solve(matX), 0);
207*bf2c3715SXin Li       VERIFY_IS_APPROX(matX, ldltup.solve(matB).eval());
208*bf2c3715SXin Li     }
209*bf2c3715SXin Li 
210*bf2c3715SXin Li     // restore
211*bf2c3715SXin Li     if(sign == -1)
212*bf2c3715SXin Li       symm = -symm;
213*bf2c3715SXin Li 
214*bf2c3715SXin Li     // check matrices coming from linear constraints with Lagrange multipliers
215*bf2c3715SXin Li     if(rows>=3)
216*bf2c3715SXin Li     {
217*bf2c3715SXin Li       SquareMatrixType A = symm;
218*bf2c3715SXin Li       Index c = internal::random<Index>(0,rows-2);
219*bf2c3715SXin Li       A.bottomRightCorner(c,c).setZero();
220*bf2c3715SXin Li       // Make sure a solution exists:
221*bf2c3715SXin Li       vecX.setRandom();
222*bf2c3715SXin Li       vecB = A * vecX;
223*bf2c3715SXin Li       vecX.setZero();
224*bf2c3715SXin Li       ldltlo.compute(A);
225*bf2c3715SXin Li       VERIFY_IS_APPROX(A, ldltlo.reconstructedMatrix());
226*bf2c3715SXin Li       vecX = ldltlo.solve(vecB);
227*bf2c3715SXin Li       VERIFY_IS_APPROX(A * vecX, vecB);
228*bf2c3715SXin Li     }
229*bf2c3715SXin Li 
230*bf2c3715SXin Li     // check non-full rank matrices
231*bf2c3715SXin Li     if(rows>=3)
232*bf2c3715SXin Li     {
233*bf2c3715SXin Li       Index r = internal::random<Index>(1,rows-1);
234*bf2c3715SXin Li       Matrix<Scalar,Dynamic,Dynamic> a = Matrix<Scalar,Dynamic,Dynamic>::Random(rows,r);
235*bf2c3715SXin Li       SquareMatrixType A = a * a.adjoint();
236*bf2c3715SXin Li       // Make sure a solution exists:
237*bf2c3715SXin Li       vecX.setRandom();
238*bf2c3715SXin Li       vecB = A * vecX;
239*bf2c3715SXin Li       vecX.setZero();
240*bf2c3715SXin Li       ldltlo.compute(A);
241*bf2c3715SXin Li       VERIFY_IS_APPROX(A, ldltlo.reconstructedMatrix());
242*bf2c3715SXin Li       vecX = ldltlo.solve(vecB);
243*bf2c3715SXin Li       VERIFY_IS_APPROX(A * vecX, vecB);
244*bf2c3715SXin Li     }
245*bf2c3715SXin Li 
246*bf2c3715SXin Li     // check matrices with a wide spectrum
247*bf2c3715SXin Li     if(rows>=3)
248*bf2c3715SXin Li     {
249*bf2c3715SXin Li       using std::pow;
250*bf2c3715SXin Li       using std::sqrt;
251*bf2c3715SXin Li       RealScalar s = (std::min)(16,std::numeric_limits<RealScalar>::max_exponent10/8);
252*bf2c3715SXin Li       Matrix<Scalar,Dynamic,Dynamic> a = Matrix<Scalar,Dynamic,Dynamic>::Random(rows,rows);
253*bf2c3715SXin Li       Matrix<RealScalar,Dynamic,1> d =  Matrix<RealScalar,Dynamic,1>::Random(rows);
254*bf2c3715SXin Li       for(Index k=0; k<rows; ++k)
255*bf2c3715SXin Li         d(k) = d(k)*pow(RealScalar(10),internal::random<RealScalar>(-s,s));
256*bf2c3715SXin Li       SquareMatrixType A = a * d.asDiagonal() * a.adjoint();
257*bf2c3715SXin Li       // Make sure a solution exists:
258*bf2c3715SXin Li       vecX.setRandom();
259*bf2c3715SXin Li       vecB = A * vecX;
260*bf2c3715SXin Li       vecX.setZero();
261*bf2c3715SXin Li       ldltlo.compute(A);
262*bf2c3715SXin Li       VERIFY_IS_APPROX(A, ldltlo.reconstructedMatrix());
263*bf2c3715SXin Li       vecX = ldltlo.solve(vecB);
264*bf2c3715SXin Li 
265*bf2c3715SXin Li       if(ldltlo.vectorD().real().cwiseAbs().minCoeff()>RealScalar(0))
266*bf2c3715SXin Li       {
267*bf2c3715SXin Li         VERIFY_IS_APPROX(A * vecX,vecB);
268*bf2c3715SXin Li       }
269*bf2c3715SXin Li       else
270*bf2c3715SXin Li       {
271*bf2c3715SXin Li         RealScalar large_tol =  sqrt(test_precision<RealScalar>());
272*bf2c3715SXin Li         VERIFY((A * vecX).isApprox(vecB, large_tol));
273*bf2c3715SXin Li 
274*bf2c3715SXin Li         ++g_test_level;
275*bf2c3715SXin Li         VERIFY_IS_APPROX(A * vecX,vecB);
276*bf2c3715SXin Li         --g_test_level;
277*bf2c3715SXin Li       }
278*bf2c3715SXin Li     }
279*bf2c3715SXin Li   }
280*bf2c3715SXin Li 
281*bf2c3715SXin Li   // update/downdate
282*bf2c3715SXin Li   CALL_SUBTEST(( test_chol_update<SquareMatrixType,LLT>(symm)  ));
283*bf2c3715SXin Li   CALL_SUBTEST(( test_chol_update<SquareMatrixType,LDLT>(symm) ));
284*bf2c3715SXin Li }
285*bf2c3715SXin Li 
cholesky_cplx(const MatrixType & m)286*bf2c3715SXin Li template<typename MatrixType> void cholesky_cplx(const MatrixType& m)
287*bf2c3715SXin Li {
288*bf2c3715SXin Li   // classic test
289*bf2c3715SXin Li   cholesky(m);
290*bf2c3715SXin Li 
291*bf2c3715SXin Li   // test mixing real/scalar types
292*bf2c3715SXin Li 
293*bf2c3715SXin Li   Index rows = m.rows();
294*bf2c3715SXin Li   Index cols = m.cols();
295*bf2c3715SXin Li 
296*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
297*bf2c3715SXin Li   typedef typename NumTraits<Scalar>::Real RealScalar;
298*bf2c3715SXin Li   typedef Matrix<RealScalar, MatrixType::RowsAtCompileTime, MatrixType::RowsAtCompileTime> RealMatrixType;
299*bf2c3715SXin Li   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
300*bf2c3715SXin Li 
301*bf2c3715SXin Li   RealMatrixType a0 = RealMatrixType::Random(rows,cols);
302*bf2c3715SXin Li   VectorType vecB = VectorType::Random(rows), vecX(rows);
303*bf2c3715SXin Li   MatrixType matB = MatrixType::Random(rows,cols), matX(rows,cols);
304*bf2c3715SXin Li   RealMatrixType symm =  a0 * a0.adjoint();
305*bf2c3715SXin Li   // let's make sure the matrix is not singular or near singular
306*bf2c3715SXin Li   for (int k=0; k<3; ++k)
307*bf2c3715SXin Li   {
308*bf2c3715SXin Li     RealMatrixType a1 = RealMatrixType::Random(rows,cols);
309*bf2c3715SXin Li     symm += a1 * a1.adjoint();
310*bf2c3715SXin Li   }
311*bf2c3715SXin Li 
312*bf2c3715SXin Li   {
313*bf2c3715SXin Li     RealMatrixType symmLo = symm.template triangularView<Lower>();
314*bf2c3715SXin Li 
315*bf2c3715SXin Li     LLT<RealMatrixType,Lower> chollo(symmLo);
316*bf2c3715SXin Li     VERIFY_IS_APPROX(symm, chollo.reconstructedMatrix());
317*bf2c3715SXin Li 
318*bf2c3715SXin Li     check_solverbase<VectorType, VectorType>(symm, chollo, rows, rows, 1);
319*bf2c3715SXin Li     //check_solverbase<MatrixType, MatrixType>(symm, chollo, rows, cols, rows);
320*bf2c3715SXin Li   }
321*bf2c3715SXin Li 
322*bf2c3715SXin Li   // LDLT
323*bf2c3715SXin Li   {
324*bf2c3715SXin Li     int sign = internal::random<int>()%2 ? 1 : -1;
325*bf2c3715SXin Li 
326*bf2c3715SXin Li     if(sign == -1)
327*bf2c3715SXin Li     {
328*bf2c3715SXin Li       symm = -symm; // test a negative matrix
329*bf2c3715SXin Li     }
330*bf2c3715SXin Li 
331*bf2c3715SXin Li     RealMatrixType symmLo = symm.template triangularView<Lower>();
332*bf2c3715SXin Li 
333*bf2c3715SXin Li     LDLT<RealMatrixType,Lower> ldltlo(symmLo);
334*bf2c3715SXin Li     VERIFY(ldltlo.info()==Success);
335*bf2c3715SXin Li     VERIFY_IS_APPROX(symm, ldltlo.reconstructedMatrix());
336*bf2c3715SXin Li 
337*bf2c3715SXin Li     check_solverbase<VectorType, VectorType>(symm, ldltlo, rows, rows, 1);
338*bf2c3715SXin Li     //check_solverbase<MatrixType, MatrixType>(symm, ldltlo, rows, cols, rows);
339*bf2c3715SXin Li   }
340*bf2c3715SXin Li }
341*bf2c3715SXin Li 
342*bf2c3715SXin Li // regression test for bug 241
cholesky_bug241(const MatrixType & m)343*bf2c3715SXin Li template<typename MatrixType> void cholesky_bug241(const MatrixType& m)
344*bf2c3715SXin Li {
345*bf2c3715SXin Li   eigen_assert(m.rows() == 2 && m.cols() == 2);
346*bf2c3715SXin Li 
347*bf2c3715SXin Li   typedef typename MatrixType::Scalar Scalar;
348*bf2c3715SXin Li   typedef Matrix<Scalar, MatrixType::RowsAtCompileTime, 1> VectorType;
349*bf2c3715SXin Li 
350*bf2c3715SXin Li   MatrixType matA;
351*bf2c3715SXin Li   matA << 1, 1, 1, 1;
352*bf2c3715SXin Li   VectorType vecB;
353*bf2c3715SXin Li   vecB << 1, 1;
354*bf2c3715SXin Li   VectorType vecX = matA.ldlt().solve(vecB);
355*bf2c3715SXin Li   VERIFY_IS_APPROX(matA * vecX, vecB);
356*bf2c3715SXin Li }
357*bf2c3715SXin Li 
358*bf2c3715SXin Li // LDLT is not guaranteed to work for indefinite matrices, but happens to work fine if matrix is diagonal.
359*bf2c3715SXin Li // This test checks that LDLT reports correctly that matrix is indefinite.
360*bf2c3715SXin Li // See http://forum.kde.org/viewtopic.php?f=74&t=106942 and bug 736
cholesky_definiteness(const MatrixType & m)361*bf2c3715SXin Li template<typename MatrixType> void cholesky_definiteness(const MatrixType& m)
362*bf2c3715SXin Li {
363*bf2c3715SXin Li   eigen_assert(m.rows() == 2 && m.cols() == 2);
364*bf2c3715SXin Li   MatrixType mat;
365*bf2c3715SXin Li   LDLT<MatrixType> ldlt(2);
366*bf2c3715SXin Li 
367*bf2c3715SXin Li   {
368*bf2c3715SXin Li     mat << 1, 0, 0, -1;
369*bf2c3715SXin Li     ldlt.compute(mat);
370*bf2c3715SXin Li     VERIFY(ldlt.info()==Success);
371*bf2c3715SXin Li     VERIFY(!ldlt.isNegative());
372*bf2c3715SXin Li     VERIFY(!ldlt.isPositive());
373*bf2c3715SXin Li     VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
374*bf2c3715SXin Li   }
375*bf2c3715SXin Li   {
376*bf2c3715SXin Li     mat << 1, 2, 2, 1;
377*bf2c3715SXin Li     ldlt.compute(mat);
378*bf2c3715SXin Li     VERIFY(ldlt.info()==Success);
379*bf2c3715SXin Li     VERIFY(!ldlt.isNegative());
380*bf2c3715SXin Li     VERIFY(!ldlt.isPositive());
381*bf2c3715SXin Li     VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
382*bf2c3715SXin Li   }
383*bf2c3715SXin Li   {
384*bf2c3715SXin Li     mat << 0, 0, 0, 0;
385*bf2c3715SXin Li     ldlt.compute(mat);
386*bf2c3715SXin Li     VERIFY(ldlt.info()==Success);
387*bf2c3715SXin Li     VERIFY(ldlt.isNegative());
388*bf2c3715SXin Li     VERIFY(ldlt.isPositive());
389*bf2c3715SXin Li     VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
390*bf2c3715SXin Li   }
391*bf2c3715SXin Li   {
392*bf2c3715SXin Li     mat << 0, 0, 0, 1;
393*bf2c3715SXin Li     ldlt.compute(mat);
394*bf2c3715SXin Li     VERIFY(ldlt.info()==Success);
395*bf2c3715SXin Li     VERIFY(!ldlt.isNegative());
396*bf2c3715SXin Li     VERIFY(ldlt.isPositive());
397*bf2c3715SXin Li     VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
398*bf2c3715SXin Li   }
399*bf2c3715SXin Li   {
400*bf2c3715SXin Li     mat << -1, 0, 0, 0;
401*bf2c3715SXin Li     ldlt.compute(mat);
402*bf2c3715SXin Li     VERIFY(ldlt.info()==Success);
403*bf2c3715SXin Li     VERIFY(ldlt.isNegative());
404*bf2c3715SXin Li     VERIFY(!ldlt.isPositive());
405*bf2c3715SXin Li     VERIFY_IS_APPROX(mat,ldlt.reconstructedMatrix());
406*bf2c3715SXin Li   }
407*bf2c3715SXin Li }
408*bf2c3715SXin Li 
409*bf2c3715SXin Li template<typename>
cholesky_faillure_cases()410*bf2c3715SXin Li void cholesky_faillure_cases()
411*bf2c3715SXin Li {
412*bf2c3715SXin Li   MatrixXd mat;
413*bf2c3715SXin Li   LDLT<MatrixXd> ldlt;
414*bf2c3715SXin Li 
415*bf2c3715SXin Li   {
416*bf2c3715SXin Li     mat.resize(2,2);
417*bf2c3715SXin Li     mat << 0, 1, 1, 0;
418*bf2c3715SXin Li     ldlt.compute(mat);
419*bf2c3715SXin Li     VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
420*bf2c3715SXin Li     VERIFY(ldlt.info()==NumericalIssue);
421*bf2c3715SXin Li   }
422*bf2c3715SXin Li #if (!EIGEN_ARCH_i386) || defined(EIGEN_VECTORIZE_SSE2)
423*bf2c3715SXin Li   {
424*bf2c3715SXin Li     mat.resize(3,3);
425*bf2c3715SXin Li     mat << -1, -3, 3,
426*bf2c3715SXin Li            -3, -8.9999999999999999999, 1,
427*bf2c3715SXin Li             3, 1, 0;
428*bf2c3715SXin Li     ldlt.compute(mat);
429*bf2c3715SXin Li     VERIFY(ldlt.info()==NumericalIssue);
430*bf2c3715SXin Li     VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
431*bf2c3715SXin Li   }
432*bf2c3715SXin Li #endif
433*bf2c3715SXin Li   {
434*bf2c3715SXin Li     mat.resize(3,3);
435*bf2c3715SXin Li     mat <<  1, 2, 3,
436*bf2c3715SXin Li             2, 4, 1,
437*bf2c3715SXin Li             3, 1, 0;
438*bf2c3715SXin Li     ldlt.compute(mat);
439*bf2c3715SXin Li     VERIFY(ldlt.info()==NumericalIssue);
440*bf2c3715SXin Li     VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
441*bf2c3715SXin Li   }
442*bf2c3715SXin Li 
443*bf2c3715SXin Li   {
444*bf2c3715SXin Li     mat.resize(8,8);
445*bf2c3715SXin Li     mat <<  0.1, 0, -0.1, 0, 0, 0, 1, 0,
446*bf2c3715SXin Li             0, 4.24667, 0, 2.00333, 0, 0, 0, 0,
447*bf2c3715SXin Li             -0.1, 0, 0.2, 0, -0.1, 0, 0, 0,
448*bf2c3715SXin Li             0, 2.00333, 0, 8.49333, 0, 2.00333, 0, 0,
449*bf2c3715SXin Li             0, 0, -0.1, 0, 0.1, 0, 0, 1,
450*bf2c3715SXin Li             0, 0, 0, 2.00333, 0, 4.24667, 0, 0,
451*bf2c3715SXin Li             1, 0, 0, 0, 0, 0, 0, 0,
452*bf2c3715SXin Li             0, 0, 0, 0, 1, 0, 0, 0;
453*bf2c3715SXin Li     ldlt.compute(mat);
454*bf2c3715SXin Li     VERIFY(ldlt.info()==NumericalIssue);
455*bf2c3715SXin Li     VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
456*bf2c3715SXin Li   }
457*bf2c3715SXin Li 
458*bf2c3715SXin Li   // bug 1479
459*bf2c3715SXin Li   {
460*bf2c3715SXin Li     mat.resize(4,4);
461*bf2c3715SXin Li     mat <<  1, 2, 0, 1,
462*bf2c3715SXin Li             2, 4, 0, 2,
463*bf2c3715SXin Li             0, 0, 0, 1,
464*bf2c3715SXin Li             1, 2, 1, 1;
465*bf2c3715SXin Li     ldlt.compute(mat);
466*bf2c3715SXin Li     VERIFY(ldlt.info()==NumericalIssue);
467*bf2c3715SXin Li     VERIFY_IS_NOT_APPROX(mat,ldlt.reconstructedMatrix());
468*bf2c3715SXin Li   }
469*bf2c3715SXin Li }
470*bf2c3715SXin Li 
cholesky_verify_assert()471*bf2c3715SXin Li template<typename MatrixType> void cholesky_verify_assert()
472*bf2c3715SXin Li {
473*bf2c3715SXin Li   MatrixType tmp;
474*bf2c3715SXin Li 
475*bf2c3715SXin Li   LLT<MatrixType> llt;
476*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(llt.matrixL())
477*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(llt.matrixU())
478*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(llt.solve(tmp))
479*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(llt.transpose().solve(tmp))
480*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(llt.adjoint().solve(tmp))
481*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(llt.solveInPlace(tmp))
482*bf2c3715SXin Li 
483*bf2c3715SXin Li   LDLT<MatrixType> ldlt;
484*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(ldlt.matrixL())
485*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(ldlt.transpositionsP())
486*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(ldlt.vectorD())
487*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(ldlt.isPositive())
488*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(ldlt.isNegative())
489*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(ldlt.solve(tmp))
490*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(ldlt.transpose().solve(tmp))
491*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(ldlt.adjoint().solve(tmp))
492*bf2c3715SXin Li   VERIFY_RAISES_ASSERT(ldlt.solveInPlace(tmp))
493*bf2c3715SXin Li }
494*bf2c3715SXin Li 
EIGEN_DECLARE_TEST(cholesky)495*bf2c3715SXin Li EIGEN_DECLARE_TEST(cholesky)
496*bf2c3715SXin Li {
497*bf2c3715SXin Li   int s = 0;
498*bf2c3715SXin Li   for(int i = 0; i < g_repeat; i++) {
499*bf2c3715SXin Li     CALL_SUBTEST_1( cholesky(Matrix<double,1,1>()) );
500*bf2c3715SXin Li     CALL_SUBTEST_3( cholesky(Matrix2d()) );
501*bf2c3715SXin Li     CALL_SUBTEST_3( cholesky_bug241(Matrix2d()) );
502*bf2c3715SXin Li     CALL_SUBTEST_3( cholesky_definiteness(Matrix2d()) );
503*bf2c3715SXin Li     CALL_SUBTEST_4( cholesky(Matrix3f()) );
504*bf2c3715SXin Li     CALL_SUBTEST_5( cholesky(Matrix4d()) );
505*bf2c3715SXin Li 
506*bf2c3715SXin Li     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE);
507*bf2c3715SXin Li     CALL_SUBTEST_2( cholesky(MatrixXd(s,s)) );
508*bf2c3715SXin Li     TEST_SET_BUT_UNUSED_VARIABLE(s)
509*bf2c3715SXin Li 
510*bf2c3715SXin Li     s = internal::random<int>(1,EIGEN_TEST_MAX_SIZE/2);
511*bf2c3715SXin Li     CALL_SUBTEST_6( cholesky_cplx(MatrixXcd(s,s)) );
512*bf2c3715SXin Li     TEST_SET_BUT_UNUSED_VARIABLE(s)
513*bf2c3715SXin Li   }
514*bf2c3715SXin Li   // empty matrix, regression test for Bug 785:
515*bf2c3715SXin Li   CALL_SUBTEST_2( cholesky(MatrixXd(0,0)) );
516*bf2c3715SXin Li 
517*bf2c3715SXin Li   // This does not work yet:
518*bf2c3715SXin Li   // CALL_SUBTEST_2( cholesky(Matrix<double,0,0>()) );
519*bf2c3715SXin Li 
520*bf2c3715SXin Li   CALL_SUBTEST_4( cholesky_verify_assert<Matrix3f>() );
521*bf2c3715SXin Li   CALL_SUBTEST_7( cholesky_verify_assert<Matrix3d>() );
522*bf2c3715SXin Li   CALL_SUBTEST_8( cholesky_verify_assert<MatrixXf>() );
523*bf2c3715SXin Li   CALL_SUBTEST_2( cholesky_verify_assert<MatrixXd>() );
524*bf2c3715SXin Li 
525*bf2c3715SXin Li   // Test problem size constructors
526*bf2c3715SXin Li   CALL_SUBTEST_9( LLT<MatrixXf>(10) );
527*bf2c3715SXin Li   CALL_SUBTEST_9( LDLT<MatrixXf>(10) );
528*bf2c3715SXin Li 
529*bf2c3715SXin Li   CALL_SUBTEST_2( cholesky_faillure_cases<void>() );
530*bf2c3715SXin Li 
531*bf2c3715SXin Li   TEST_SET_BUT_UNUSED_VARIABLE(nb_temporaries)
532*bf2c3715SXin Li }
533