aboutsummaryrefslogtreecommitdiffhomepage
path: root/unsupported/test/sparse_lu.cpp
blob: bb0dccdb540ecc597685fe9f3c8ebe1684987484 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2008-2010 Gael Guennebaud <g.gael@free.fr>
//
// Eigen is free software; you can redistribute it and/or
// modify it under the terms of the GNU Lesser General Public
// License as published by the Free Software Foundation; either
// version 3 of the License, or (at your option) any later version.
//
// Alternatively, you can redistribute it and/or
// modify it under the terms of the GNU General Public License as
// published by the Free Software Foundation; either version 2 of
// the License, or (at your option) any later version.
//
// Eigen is distributed in the hope that it will be useful, but WITHOUT ANY
// WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
// FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License or the
// GNU General Public License for more details.
//
// You should have received a copy of the GNU Lesser General Public
// License and a copy of the GNU General Public License along with
// Eigen. If not, see <http://www.gnu.org/licenses/>.

#include "sparse.h"
#include <Eigen/SparseExtra>

#ifdef EIGEN_UMFPACK_SUPPORT
#include <Eigen/UmfPackSupport>
#endif

#ifdef EIGEN_SUPERLU_SUPPORT
#include <Eigen/SuperLUSupport>
#endif

template<typename Scalar> void sparse_lu(int rows, int cols)
{
  double density = std::max(8./(rows*cols), 0.01);
  typedef Matrix<Scalar,Dynamic,Dynamic> DenseMatrix;
  typedef Matrix<Scalar,Dynamic,1> DenseVector;

  DenseVector vec1 = DenseVector::Random(rows);

  std::vector<Vector2i> zeroCoords;
  std::vector<Vector2i> nonzeroCoords;

    SparseMatrix<Scalar> m2(rows, cols);
    DenseMatrix refMat2(rows, cols);

    DenseVector b = DenseVector::Random(cols);
    DenseVector refX(cols), x(cols);

    initSparse<Scalar>(density, refMat2, m2, ForceNonZeroDiag, &zeroCoords, &nonzeroCoords);

    FullPivLU<DenseMatrix> refLu(refMat2);
    refX = refLu.solve(b);
    #if defined(EIGEN_SUPERLU_SUPPORT) || defined(EIGEN_UMFPACK_SUPPORT)
    Scalar refDet = refLu.determinant();
    #endif
    x.setZero();
    // // SparseLU<SparseMatrix<Scalar> > (m2).solve(b,&x);
    // // VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LU: default");
    
    #ifdef EIGEN_UMFPACK_SUPPORT
    {
      // check solve
      x.setZero();
      SparseLU<SparseMatrix<Scalar>,UmfPack> lu(m2);
      VERIFY(lu.succeeded() && "umfpack LU decomposition failed");
      VERIFY(lu.solve(b,&x) && "umfpack LU solving failed");
      VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LU: umfpack");
      VERIFY_IS_APPROX(refDet,lu.determinant());
        // TODO check the extracted data
        //std::cerr << slu.matrixL() << "\n";
    }
    #endif
    
    #ifdef EIGEN_SUPERLU_SUPPORT
    // legacy, deprecated API
    {
      x.setZero();
      SparseLU<SparseMatrix<Scalar>,SuperLULegacy> slu(m2);
      if (slu.succeeded())
      {
        DenseVector oldb = b;
        if (slu.solve(b,&x)) {
          VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "LU: SuperLU");
        }
        else
          std::cerr << "super lu solving failed\n";
        VERIFY(oldb.isApprox(b) && "the rhs should not be modified!");
        
        // std::cerr << refDet << " == " << slu.determinant() << "\n";
        if (slu.solve(b, &x, SvTranspose)) {
          VERIFY(b.isApprox(m2.transpose() * x, test_precision<Scalar>()));
        }
        else
          std::cerr << "super lu solving failed\n";

        if (slu.solve(b, &x, SvAdjoint)) {
         VERIFY(b.isApprox(m2.adjoint() * x, test_precision<Scalar>()));
        }
        else
          std::cerr << "super lu solving failed\n";

        if (!NumTraits<Scalar>::IsComplex) {
          VERIFY_IS_APPROX(refDet,slu.determinant()); // FIXME det is not very stable for complex
        }
      }
      else
        std::cerr << "super lu factorize failed\n";
    }
    
    // New API
    {
      x.setZero();
      SuperLU<SparseMatrix<Scalar> > slu(m2);
      if (slu.info() == Success)
      {
        DenseVector oldb = b;
        x = slu.solve(b);
        VERIFY(oldb.isApprox(b) && "the rhs should not be modified!");
        if (slu.info() == Success) {
          VERIFY(refX.isApprox(x,test_precision<Scalar>()) && "SuperLU");
        }
        else
          std::cerr << "super lu solving failed\n";

        if (!NumTraits<Scalar>::IsComplex) {
          VERIFY_IS_APPROX(refDet,slu.determinant()); // FIXME det is not very stable for complex
        }
      }
      else
        std::cerr << "super lu factorize failed\n";
    }
    #endif
    
}

void test_sparse_lu()
{
  for(int i = 0; i < g_repeat; i++) {
    CALL_SUBTEST_1(sparse_lu<double>(8, 8) );
    int s = internal::random<int>(1,300);
    CALL_SUBTEST_2(sparse_lu<std::complex<double> >(s,s) );
    CALL_SUBTEST_1(sparse_lu<double>(s,s) );
  }
}