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) 2015-2016 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 // #define EIGEN_DONT_VECTORIZE
10*bf2c3715SXin Li // #define EIGEN_MAX_ALIGN_BYTES 0
11*bf2c3715SXin Li #include "sparse_solver.h"
12*bf2c3715SXin Li #include <Eigen/IterativeLinearSolvers>
13*bf2c3715SXin Li #include <unsupported/Eigen/IterativeSolvers>
14*bf2c3715SXin Li
test_incomplete_cholesky_T()15*bf2c3715SXin Li template<typename T, typename I_> void test_incomplete_cholesky_T()
16*bf2c3715SXin Li {
17*bf2c3715SXin Li typedef SparseMatrix<T,0,I_> SparseMatrixType;
18*bf2c3715SXin Li ConjugateGradient<SparseMatrixType, Lower, IncompleteCholesky<T, Lower, AMDOrdering<I_> > > cg_illt_lower_amd;
19*bf2c3715SXin Li ConjugateGradient<SparseMatrixType, Lower, IncompleteCholesky<T, Lower, NaturalOrdering<I_> > > cg_illt_lower_nat;
20*bf2c3715SXin Li ConjugateGradient<SparseMatrixType, Upper, IncompleteCholesky<T, Upper, AMDOrdering<I_> > > cg_illt_upper_amd;
21*bf2c3715SXin Li ConjugateGradient<SparseMatrixType, Upper, IncompleteCholesky<T, Upper, NaturalOrdering<I_> > > cg_illt_upper_nat;
22*bf2c3715SXin Li ConjugateGradient<SparseMatrixType, Upper|Lower, IncompleteCholesky<T, Lower, AMDOrdering<I_> > > cg_illt_uplo_amd;
23*bf2c3715SXin Li
24*bf2c3715SXin Li
25*bf2c3715SXin Li CALL_SUBTEST( check_sparse_spd_solving(cg_illt_lower_amd) );
26*bf2c3715SXin Li CALL_SUBTEST( check_sparse_spd_solving(cg_illt_lower_nat) );
27*bf2c3715SXin Li CALL_SUBTEST( check_sparse_spd_solving(cg_illt_upper_amd) );
28*bf2c3715SXin Li CALL_SUBTEST( check_sparse_spd_solving(cg_illt_upper_nat) );
29*bf2c3715SXin Li CALL_SUBTEST( check_sparse_spd_solving(cg_illt_uplo_amd) );
30*bf2c3715SXin Li }
31*bf2c3715SXin Li
32*bf2c3715SXin Li template<int>
bug1150()33*bf2c3715SXin Li void bug1150()
34*bf2c3715SXin Li {
35*bf2c3715SXin Li // regression for bug 1150
36*bf2c3715SXin Li for(int N = 1; N<20; ++N)
37*bf2c3715SXin Li {
38*bf2c3715SXin Li Eigen::MatrixXd b( N, N );
39*bf2c3715SXin Li b.setOnes();
40*bf2c3715SXin Li
41*bf2c3715SXin Li Eigen::SparseMatrix<double> m( N, N );
42*bf2c3715SXin Li m.reserve(Eigen::VectorXi::Constant(N,4));
43*bf2c3715SXin Li for( int i = 0; i < N; ++i )
44*bf2c3715SXin Li {
45*bf2c3715SXin Li m.insert( i, i ) = 1;
46*bf2c3715SXin Li m.coeffRef( i, i / 2 ) = 2;
47*bf2c3715SXin Li m.coeffRef( i, i / 3 ) = 2;
48*bf2c3715SXin Li m.coeffRef( i, i / 4 ) = 2;
49*bf2c3715SXin Li }
50*bf2c3715SXin Li
51*bf2c3715SXin Li Eigen::SparseMatrix<double> A;
52*bf2c3715SXin Li A = m * m.transpose();
53*bf2c3715SXin Li
54*bf2c3715SXin Li Eigen::ConjugateGradient<Eigen::SparseMatrix<double>,
55*bf2c3715SXin Li Eigen::Lower | Eigen::Upper,
56*bf2c3715SXin Li Eigen::IncompleteCholesky<double> > solver( A );
57*bf2c3715SXin Li VERIFY(solver.preconditioner().info() == Eigen::Success);
58*bf2c3715SXin Li VERIFY(solver.info() == Eigen::Success);
59*bf2c3715SXin Li }
60*bf2c3715SXin Li }
61*bf2c3715SXin Li
EIGEN_DECLARE_TEST(incomplete_cholesky)62*bf2c3715SXin Li EIGEN_DECLARE_TEST(incomplete_cholesky)
63*bf2c3715SXin Li {
64*bf2c3715SXin Li CALL_SUBTEST_1(( test_incomplete_cholesky_T<double,int>() ));
65*bf2c3715SXin Li CALL_SUBTEST_2(( test_incomplete_cholesky_T<std::complex<double>, int>() ));
66*bf2c3715SXin Li CALL_SUBTEST_3(( test_incomplete_cholesky_T<double,long int>() ));
67*bf2c3715SXin Li
68*bf2c3715SXin Li CALL_SUBTEST_1(( bug1150<0>() ));
69*bf2c3715SXin Li }
70