1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2012 Désiré Nuentsa-Wakam <desire.nuentsa_wakam@inria.fr>
5// Copyright (C) 2012 Gael Guennebaud <gael.guennebaud@inria.fr>
6//
7// This Source Code Form is subject to the terms of the Mozilla
8// Public License v. 2.0. If a copy of the MPL was not distributed
9// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
10
11
12#ifndef EIGEN_SPARSE_LU_H
13#define EIGEN_SPARSE_LU_H
14
15namespace Eigen {
16
17template <typename _MatrixType, typename _OrderingType = COLAMDOrdering<typename _MatrixType::Index> > class SparseLU;
18template <typename MappedSparseMatrixType> struct SparseLUMatrixLReturnType;
19template <typename MatrixLType, typename MatrixUType> struct SparseLUMatrixUReturnType;
20
21/** \ingroup SparseLU_Module
22  * \class SparseLU
23  *
24  * \brief Sparse supernodal LU factorization for general matrices
25  *
26  * This class implements the supernodal LU factorization for general matrices.
27  * It uses the main techniques from the sequential SuperLU package
28  * (http://crd-legacy.lbl.gov/~xiaoye/SuperLU/). It handles transparently real
29  * and complex arithmetics with single and double precision, depending on the
30  * scalar type of your input matrix.
31  * The code has been optimized to provide BLAS-3 operations during supernode-panel updates.
32  * It benefits directly from the built-in high-performant Eigen BLAS routines.
33  * Moreover, when the size of a supernode is very small, the BLAS calls are avoided to
34  * enable a better optimization from the compiler. For best performance,
35  * you should compile it with NDEBUG flag to avoid the numerous bounds checking on vectors.
36  *
37  * An important parameter of this class is the ordering method. It is used to reorder the columns
38  * (and eventually the rows) of the matrix to reduce the number of new elements that are created during
39  * numerical factorization. The cheapest method available is COLAMD.
40  * See  \link OrderingMethods_Module the OrderingMethods module \endlink for the list of
41  * built-in and external ordering methods.
42  *
43  * Simple example with key steps
44  * \code
45  * VectorXd x(n), b(n);
46  * SparseMatrix<double, ColMajor> A;
47  * SparseLU<SparseMatrix<scalar, ColMajor>, COLAMDOrdering<Index> >   solver;
48  * // fill A and b;
49  * // Compute the ordering permutation vector from the structural pattern of A
50  * solver.analyzePattern(A);
51  * // Compute the numerical factorization
52  * solver.factorize(A);
53  * //Use the factors to solve the linear system
54  * x = solver.solve(b);
55  * \endcode
56  *
57  * \warning The input matrix A should be in a \b compressed and \b column-major form.
58  * Otherwise an expensive copy will be made. You can call the inexpensive makeCompressed() to get a compressed matrix.
59  *
60  * \note Unlike the initial SuperLU implementation, there is no step to equilibrate the matrix.
61  * For badly scaled matrices, this step can be useful to reduce the pivoting during factorization.
62  * If this is the case for your matrices, you can try the basic scaling method at
63  *  "unsupported/Eigen/src/IterativeSolvers/Scaling.h"
64  *
65  * \tparam _MatrixType The type of the sparse matrix. It must be a column-major SparseMatrix<>
66  * \tparam _OrderingType The ordering method to use, either AMD, COLAMD or METIS. Default is COLMAD
67  *
68  *
69  * \sa \ref TutorialSparseDirectSolvers
70  * \sa \ref OrderingMethods_Module
71  */
72template <typename _MatrixType, typename _OrderingType>
73class SparseLU : public internal::SparseLUImpl<typename _MatrixType::Scalar, typename _MatrixType::Index>
74{
75  public:
76    typedef _MatrixType MatrixType;
77    typedef _OrderingType OrderingType;
78    typedef typename MatrixType::Scalar Scalar;
79    typedef typename MatrixType::RealScalar RealScalar;
80    typedef typename MatrixType::Index Index;
81    typedef SparseMatrix<Scalar,ColMajor,Index> NCMatrix;
82    typedef internal::MappedSuperNodalMatrix<Scalar, Index> SCMatrix;
83    typedef Matrix<Scalar,Dynamic,1> ScalarVector;
84    typedef Matrix<Index,Dynamic,1> IndexVector;
85    typedef PermutationMatrix<Dynamic, Dynamic, Index> PermutationType;
86    typedef internal::SparseLUImpl<Scalar, Index> Base;
87
88  public:
89    SparseLU():m_isInitialized(true),m_lastError(""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(false),m_diagpivotthresh(1.0),m_detPermR(1)
90    {
91      initperfvalues();
92    }
93    SparseLU(const MatrixType& matrix):m_isInitialized(true),m_lastError(""),m_Ustore(0,0,0,0,0,0),m_symmetricmode(false),m_diagpivotthresh(1.0),m_detPermR(1)
94    {
95      initperfvalues();
96      compute(matrix);
97    }
98
99    ~SparseLU()
100    {
101      // Free all explicit dynamic pointers
102    }
103
104    void analyzePattern (const MatrixType& matrix);
105    void factorize (const MatrixType& matrix);
106    void simplicialfactorize(const MatrixType& matrix);
107
108    /**
109      * Compute the symbolic and numeric factorization of the input sparse matrix.
110      * The input matrix should be in column-major storage.
111      */
112    void compute (const MatrixType& matrix)
113    {
114      // Analyze
115      analyzePattern(matrix);
116      //Factorize
117      factorize(matrix);
118    }
119
120    inline Index rows() const { return m_mat.rows(); }
121    inline Index cols() const { return m_mat.cols(); }
122    /** Indicate that the pattern of the input matrix is symmetric */
123    void isSymmetric(bool sym)
124    {
125      m_symmetricmode = sym;
126    }
127
128    /** \returns an expression of the matrix L, internally stored as supernodes
129      * The only operation available with this expression is the triangular solve
130      * \code
131      * y = b; matrixL().solveInPlace(y);
132      * \endcode
133      */
134    SparseLUMatrixLReturnType<SCMatrix> matrixL() const
135    {
136      return SparseLUMatrixLReturnType<SCMatrix>(m_Lstore);
137    }
138    /** \returns an expression of the matrix U,
139      * The only operation available with this expression is the triangular solve
140      * \code
141      * y = b; matrixU().solveInPlace(y);
142      * \endcode
143      */
144    SparseLUMatrixUReturnType<SCMatrix,MappedSparseMatrix<Scalar,ColMajor,Index> > matrixU() const
145    {
146      return SparseLUMatrixUReturnType<SCMatrix, MappedSparseMatrix<Scalar,ColMajor,Index> >(m_Lstore, m_Ustore);
147    }
148
149    /**
150      * \returns a reference to the row matrix permutation \f$ P_r \f$ such that \f$P_r A P_c^T = L U\f$
151      * \sa colsPermutation()
152      */
153    inline const PermutationType& rowsPermutation() const
154    {
155      return m_perm_r;
156    }
157    /**
158      * \returns a reference to the column matrix permutation\f$ P_c^T \f$ such that \f$P_r A P_c^T = L U\f$
159      * \sa rowsPermutation()
160      */
161    inline const PermutationType& colsPermutation() const
162    {
163      return m_perm_c;
164    }
165    /** Set the threshold used for a diagonal entry to be an acceptable pivot. */
166    void setPivotThreshold(const RealScalar& thresh)
167    {
168      m_diagpivotthresh = thresh;
169    }
170
171    /** \returns the solution X of \f$ A X = B \f$ using the current decomposition of A.
172      *
173      * \warning the destination matrix X in X = this->solve(B) must be colmun-major.
174      *
175      * \sa compute()
176      */
177    template<typename Rhs>
178    inline const internal::solve_retval<SparseLU, Rhs> solve(const MatrixBase<Rhs>& B) const
179    {
180      eigen_assert(m_factorizationIsOk && "SparseLU is not initialized.");
181      eigen_assert(rows()==B.rows()
182                    && "SparseLU::solve(): invalid number of rows of the right hand side matrix B");
183          return internal::solve_retval<SparseLU, Rhs>(*this, B.derived());
184    }
185
186    /** \returns the solution X of \f$ A X = B \f$ using the current decomposition of A.
187      *
188      * \sa compute()
189      */
190    template<typename Rhs>
191    inline const internal::sparse_solve_retval<SparseLU, Rhs> solve(const SparseMatrixBase<Rhs>& B) const
192    {
193      eigen_assert(m_factorizationIsOk && "SparseLU is not initialized.");
194      eigen_assert(rows()==B.rows()
195                    && "SparseLU::solve(): invalid number of rows of the right hand side matrix B");
196          return internal::sparse_solve_retval<SparseLU, Rhs>(*this, B.derived());
197    }
198
199    /** \brief Reports whether previous computation was successful.
200      *
201      * \returns \c Success if computation was succesful,
202      *          \c NumericalIssue if the LU factorization reports a problem, zero diagonal for instance
203      *          \c InvalidInput if the input matrix is invalid
204      *
205      * \sa iparm()
206      */
207    ComputationInfo info() const
208    {
209      eigen_assert(m_isInitialized && "Decomposition is not initialized.");
210      return m_info;
211    }
212
213    /**
214      * \returns A string describing the type of error
215      */
216    std::string lastErrorMessage() const
217    {
218      return m_lastError;
219    }
220
221    template<typename Rhs, typename Dest>
222    bool _solve(const MatrixBase<Rhs> &B, MatrixBase<Dest> &X_base) const
223    {
224      Dest& X(X_base.derived());
225      eigen_assert(m_factorizationIsOk && "The matrix should be factorized first");
226      EIGEN_STATIC_ASSERT((Dest::Flags&RowMajorBit)==0,
227                        THIS_METHOD_IS_ONLY_FOR_COLUMN_MAJOR_MATRICES);
228
229      // Permute the right hand side to form X = Pr*B
230      // on return, X is overwritten by the computed solution
231      X.resize(B.rows(),B.cols());
232
233      // this ugly const_cast_derived() helps to detect aliasing when applying the permutations
234      for(Index j = 0; j < B.cols(); ++j)
235        X.col(j) = rowsPermutation() * B.const_cast_derived().col(j);
236
237      //Forward substitution with L
238      this->matrixL().solveInPlace(X);
239      this->matrixU().solveInPlace(X);
240
241      // Permute back the solution
242      for (Index j = 0; j < B.cols(); ++j)
243        X.col(j) = colsPermutation().inverse() * X.col(j);
244
245      return true;
246    }
247
248    /**
249      * \returns the absolute value of the determinant of the matrix of which
250      * *this is the QR decomposition.
251      *
252      * \warning a determinant can be very big or small, so for matrices
253      * of large enough dimension, there is a risk of overflow/underflow.
254      * One way to work around that is to use logAbsDeterminant() instead.
255      *
256      * \sa logAbsDeterminant(), signDeterminant()
257      */
258     Scalar absDeterminant()
259    {
260      eigen_assert(m_factorizationIsOk && "The matrix should be factorized first.");
261      // Initialize with the determinant of the row matrix
262      Scalar det = Scalar(1.);
263      // Note that the diagonal blocks of U are stored in supernodes,
264      // which are available in the  L part :)
265      for (Index j = 0; j < this->cols(); ++j)
266      {
267        for (typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
268        {
269          if(it.index() == j)
270          {
271            det *= (std::abs)(it.value());
272            break;
273          }
274        }
275       }
276       return det;
277     }
278
279     /** \returns the natural log of the absolute value of the determinant of the matrix
280       * of which **this is the QR decomposition
281       *
282       * \note This method is useful to work around the risk of overflow/underflow that's
283       * inherent to the determinant computation.
284       *
285       * \sa absDeterminant(), signDeterminant()
286       */
287     Scalar logAbsDeterminant() const
288     {
289       eigen_assert(m_factorizationIsOk && "The matrix should be factorized first.");
290       Scalar det = Scalar(0.);
291       for (Index j = 0; j < this->cols(); ++j)
292       {
293         for (typename SCMatrix::InnerIterator it(m_Lstore, j); it; ++it)
294         {
295           if(it.row() < j) continue;
296           if(it.row() == j)
297           {
298             det += (std::log)((std::abs)(it.value()));
299             break;
300           }
301         }
302       }
303       return det;
304     }
305
306     /** \returns A number representing the sign of the determinant
307       *
308       * \sa absDeterminant(), logAbsDeterminant()
309       */
310     Scalar signDeterminant()
311     {
312       eigen_assert(m_factorizationIsOk && "The matrix should be factorized first.");
313       return Scalar(m_detPermR);
314     }
315
316  protected:
317    // Functions
318    void initperfvalues()
319    {
320      m_perfv.panel_size = 1;
321      m_perfv.relax = 1;
322      m_perfv.maxsuper = 128;
323      m_perfv.rowblk = 16;
324      m_perfv.colblk = 8;
325      m_perfv.fillfactor = 20;
326    }
327
328    // Variables
329    mutable ComputationInfo m_info;
330    bool m_isInitialized;
331    bool m_factorizationIsOk;
332    bool m_analysisIsOk;
333    std::string m_lastError;
334    NCMatrix m_mat; // The input (permuted ) matrix
335    SCMatrix m_Lstore; // The lower triangular matrix (supernodal)
336    MappedSparseMatrix<Scalar,ColMajor,Index> m_Ustore; // The upper triangular matrix
337    PermutationType m_perm_c; // Column permutation
338    PermutationType m_perm_r ; // Row permutation
339    IndexVector m_etree; // Column elimination tree
340
341    typename Base::GlobalLU_t m_glu;
342
343    // SparseLU options
344    bool m_symmetricmode;
345    // values for performance
346    internal::perfvalues<Index> m_perfv;
347    RealScalar m_diagpivotthresh; // Specifies the threshold used for a diagonal entry to be an acceptable pivot
348    Index m_nnzL, m_nnzU; // Nonzeros in L and U factors
349    Index m_detPermR; // Determinant of the coefficient matrix
350  private:
351    // Disable copy constructor
352    SparseLU (const SparseLU& );
353
354}; // End class SparseLU
355
356
357
358// Functions needed by the anaysis phase
359/**
360  * Compute the column permutation to minimize the fill-in
361  *
362  *  - Apply this permutation to the input matrix -
363  *
364  *  - Compute the column elimination tree on the permuted matrix
365  *
366  *  - Postorder the elimination tree and the column permutation
367  *
368  */
369template <typename MatrixType, typename OrderingType>
370void SparseLU<MatrixType, OrderingType>::analyzePattern(const MatrixType& mat)
371{
372
373  //TODO  It is possible as in SuperLU to compute row and columns scaling vectors to equilibrate the matrix mat.
374
375  OrderingType ord;
376  ord(mat,m_perm_c);
377
378  // Apply the permutation to the column of the input  matrix
379  //First copy the whole input matrix.
380  m_mat = mat;
381  if (m_perm_c.size()) {
382    m_mat.uncompress(); //NOTE: The effect of this command is only to create the InnerNonzeros pointers. FIXME : This vector is filled but not subsequently used.
383    //Then, permute only the column pointers
384    const Index * outerIndexPtr;
385    if (mat.isCompressed()) outerIndexPtr = mat.outerIndexPtr();
386    else
387    {
388      Index *outerIndexPtr_t = new Index[mat.cols()+1];
389      for(Index i = 0; i <= mat.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i];
390      outerIndexPtr = outerIndexPtr_t;
391    }
392    for (Index i = 0; i < mat.cols(); i++)
393    {
394      m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
395      m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
396    }
397    if(!mat.isCompressed()) delete[] outerIndexPtr;
398  }
399  // Compute the column elimination tree of the permuted matrix
400  IndexVector firstRowElt;
401  internal::coletree(m_mat, m_etree,firstRowElt);
402
403  // In symmetric mode, do not do postorder here
404  if (!m_symmetricmode) {
405    IndexVector post, iwork;
406    // Post order etree
407    internal::treePostorder(m_mat.cols(), m_etree, post);
408
409
410    // Renumber etree in postorder
411    Index m = m_mat.cols();
412    iwork.resize(m+1);
413    for (Index i = 0; i < m; ++i) iwork(post(i)) = post(m_etree(i));
414    m_etree = iwork;
415
416    // Postmultiply A*Pc by post, i.e reorder the matrix according to the postorder of the etree
417    PermutationType post_perm(m);
418    for (Index i = 0; i < m; i++)
419      post_perm.indices()(i) = post(i);
420
421    // Combine the two permutations : postorder the permutation for future use
422    if(m_perm_c.size()) {
423      m_perm_c = post_perm * m_perm_c;
424    }
425
426  } // end postordering
427
428  m_analysisIsOk = true;
429}
430
431// Functions needed by the numerical factorization phase
432
433
434/**
435  *  - Numerical factorization
436  *  - Interleaved with the symbolic factorization
437  * On exit,  info is
438  *
439  *    = 0: successful factorization
440  *
441  *    > 0: if info = i, and i is
442  *
443  *       <= A->ncol: U(i,i) is exactly zero. The factorization has
444  *          been completed, but the factor U is exactly singular,
445  *          and division by zero will occur if it is used to solve a
446  *          system of equations.
447  *
448  *       > A->ncol: number of bytes allocated when memory allocation
449  *         failure occurred, plus A->ncol. If lwork = -1, it is
450  *         the estimated amount of space needed, plus A->ncol.
451  */
452template <typename MatrixType, typename OrderingType>
453void SparseLU<MatrixType, OrderingType>::factorize(const MatrixType& matrix)
454{
455  using internal::emptyIdxLU;
456  eigen_assert(m_analysisIsOk && "analyzePattern() should be called first");
457  eigen_assert((matrix.rows() == matrix.cols()) && "Only for squared matrices");
458
459  typedef typename IndexVector::Scalar Index;
460
461
462  // Apply the column permutation computed in analyzepattern()
463  //   m_mat = matrix * m_perm_c.inverse();
464  m_mat = matrix;
465  if (m_perm_c.size())
466  {
467    m_mat.uncompress(); //NOTE: The effect of this command is only to create the InnerNonzeros pointers.
468    //Then, permute only the column pointers
469    const Index * outerIndexPtr;
470    if (matrix.isCompressed()) outerIndexPtr = matrix.outerIndexPtr();
471    else
472    {
473      Index* outerIndexPtr_t = new Index[matrix.cols()+1];
474      for(Index i = 0; i <= matrix.cols(); i++) outerIndexPtr_t[i] = m_mat.outerIndexPtr()[i];
475      outerIndexPtr = outerIndexPtr_t;
476    }
477    for (Index i = 0; i < matrix.cols(); i++)
478    {
479      m_mat.outerIndexPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i];
480      m_mat.innerNonZeroPtr()[m_perm_c.indices()(i)] = outerIndexPtr[i+1] - outerIndexPtr[i];
481    }
482    if(!matrix.isCompressed()) delete[] outerIndexPtr;
483  }
484  else
485  { //FIXME This should not be needed if the empty permutation is handled transparently
486    m_perm_c.resize(matrix.cols());
487    for(Index i = 0; i < matrix.cols(); ++i) m_perm_c.indices()(i) = i;
488  }
489
490  Index m = m_mat.rows();
491  Index n = m_mat.cols();
492  Index nnz = m_mat.nonZeros();
493  Index maxpanel = m_perfv.panel_size * m;
494  // Allocate working storage common to the factor routines
495  Index lwork = 0;
496  Index info = Base::memInit(m, n, nnz, lwork, m_perfv.fillfactor, m_perfv.panel_size, m_glu);
497  if (info)
498  {
499    m_lastError = "UNABLE TO ALLOCATE WORKING MEMORY\n\n" ;
500    m_factorizationIsOk = false;
501    return ;
502  }
503
504  // Set up pointers for integer working arrays
505  IndexVector segrep(m); segrep.setZero();
506  IndexVector parent(m); parent.setZero();
507  IndexVector xplore(m); xplore.setZero();
508  IndexVector repfnz(maxpanel);
509  IndexVector panel_lsub(maxpanel);
510  IndexVector xprune(n); xprune.setZero();
511  IndexVector marker(m*internal::LUNoMarker); marker.setZero();
512
513  repfnz.setConstant(-1);
514  panel_lsub.setConstant(-1);
515
516  // Set up pointers for scalar working arrays
517  ScalarVector dense;
518  dense.setZero(maxpanel);
519  ScalarVector tempv;
520  tempv.setZero(internal::LUnumTempV(m, m_perfv.panel_size, m_perfv.maxsuper, /*m_perfv.rowblk*/m) );
521
522  // Compute the inverse of perm_c
523  PermutationType iperm_c(m_perm_c.inverse());
524
525  // Identify initial relaxed snodes
526  IndexVector relax_end(n);
527  if ( m_symmetricmode == true )
528    Base::heap_relax_snode(n, m_etree, m_perfv.relax, marker, relax_end);
529  else
530    Base::relax_snode(n, m_etree, m_perfv.relax, marker, relax_end);
531
532
533  m_perm_r.resize(m);
534  m_perm_r.indices().setConstant(-1);
535  marker.setConstant(-1);
536  m_detPermR = 1; // Record the determinant of the row permutation
537
538  m_glu.supno(0) = emptyIdxLU; m_glu.xsup.setConstant(0);
539  m_glu.xsup(0) = m_glu.xlsub(0) = m_glu.xusub(0) = m_glu.xlusup(0) = Index(0);
540
541  // Work on one 'panel' at a time. A panel is one of the following :
542  //  (a) a relaxed supernode at the bottom of the etree, or
543  //  (b) panel_size contiguous columns, <panel_size> defined by the user
544  Index jcol;
545  IndexVector panel_histo(n);
546  Index pivrow; // Pivotal row number in the original row matrix
547  Index nseg1; // Number of segments in U-column above panel row jcol
548  Index nseg; // Number of segments in each U-column
549  Index irep;
550  Index i, k, jj;
551  for (jcol = 0; jcol < n; )
552  {
553    // Adjust panel size so that a panel won't overlap with the next relaxed snode.
554    Index panel_size = m_perfv.panel_size; // upper bound on panel width
555    for (k = jcol + 1; k < (std::min)(jcol+panel_size, n); k++)
556    {
557      if (relax_end(k) != emptyIdxLU)
558      {
559        panel_size = k - jcol;
560        break;
561      }
562    }
563    if (k == n)
564      panel_size = n - jcol;
565
566    // Symbolic outer factorization on a panel of columns
567    Base::panel_dfs(m, panel_size, jcol, m_mat, m_perm_r.indices(), nseg1, dense, panel_lsub, segrep, repfnz, xprune, marker, parent, xplore, m_glu);
568
569    // Numeric sup-panel updates in topological order
570    Base::panel_bmod(m, panel_size, jcol, nseg1, dense, tempv, segrep, repfnz, m_glu);
571
572    // Sparse LU within the panel, and below the panel diagonal
573    for ( jj = jcol; jj< jcol + panel_size; jj++)
574    {
575      k = (jj - jcol) * m; // Column index for w-wide arrays
576
577      nseg = nseg1; // begin after all the panel segments
578      //Depth-first-search for the current column
579      VectorBlock<IndexVector> panel_lsubk(panel_lsub, k, m);
580      VectorBlock<IndexVector> repfnz_k(repfnz, k, m);
581      info = Base::column_dfs(m, jj, m_perm_r.indices(), m_perfv.maxsuper, nseg, panel_lsubk, segrep, repfnz_k, xprune, marker, parent, xplore, m_glu);
582      if ( info )
583      {
584        m_lastError =  "UNABLE TO EXPAND MEMORY IN COLUMN_DFS() ";
585        m_info = NumericalIssue;
586        m_factorizationIsOk = false;
587        return;
588      }
589      // Numeric updates to this column
590      VectorBlock<ScalarVector> dense_k(dense, k, m);
591      VectorBlock<IndexVector> segrep_k(segrep, nseg1, m-nseg1);
592      info = Base::column_bmod(jj, (nseg - nseg1), dense_k, tempv, segrep_k, repfnz_k, jcol, m_glu);
593      if ( info )
594      {
595        m_lastError = "UNABLE TO EXPAND MEMORY IN COLUMN_BMOD() ";
596        m_info = NumericalIssue;
597        m_factorizationIsOk = false;
598        return;
599      }
600
601      // Copy the U-segments to ucol(*)
602      info = Base::copy_to_ucol(jj, nseg, segrep, repfnz_k ,m_perm_r.indices(), dense_k, m_glu);
603      if ( info )
604      {
605        m_lastError = "UNABLE TO EXPAND MEMORY IN COPY_TO_UCOL() ";
606        m_info = NumericalIssue;
607        m_factorizationIsOk = false;
608        return;
609      }
610
611      // Form the L-segment
612      info = Base::pivotL(jj, m_diagpivotthresh, m_perm_r.indices(), iperm_c.indices(), pivrow, m_glu);
613      if ( info )
614      {
615        m_lastError = "THE MATRIX IS STRUCTURALLY SINGULAR ... ZERO COLUMN AT ";
616        std::ostringstream returnInfo;
617        returnInfo << info;
618        m_lastError += returnInfo.str();
619        m_info = NumericalIssue;
620        m_factorizationIsOk = false;
621        return;
622      }
623
624      // Update the determinant of the row permutation matrix
625      if (pivrow != jj) m_detPermR *= -1;
626
627      // Prune columns (0:jj-1) using column jj
628      Base::pruneL(jj, m_perm_r.indices(), pivrow, nseg, segrep, repfnz_k, xprune, m_glu);
629
630      // Reset repfnz for this column
631      for (i = 0; i < nseg; i++)
632      {
633        irep = segrep(i);
634        repfnz_k(irep) = emptyIdxLU;
635      }
636    } // end SparseLU within the panel
637    jcol += panel_size;  // Move to the next panel
638  } // end for -- end elimination
639
640  // Count the number of nonzeros in factors
641  Base::countnz(n, m_nnzL, m_nnzU, m_glu);
642  // Apply permutation  to the L subscripts
643  Base::fixupL(n, m_perm_r.indices(), m_glu);
644
645  // Create supernode matrix L
646  m_Lstore.setInfos(m, n, m_glu.lusup, m_glu.xlusup, m_glu.lsub, m_glu.xlsub, m_glu.supno, m_glu.xsup);
647  // Create the column major upper sparse matrix  U;
648  new (&m_Ustore) MappedSparseMatrix<Scalar, ColMajor, Index> ( m, n, m_nnzU, m_glu.xusub.data(), m_glu.usub.data(), m_glu.ucol.data() );
649
650  m_info = Success;
651  m_factorizationIsOk = true;
652}
653
654template<typename MappedSupernodalType>
655struct SparseLUMatrixLReturnType : internal::no_assignment_operator
656{
657  typedef typename MappedSupernodalType::Index Index;
658  typedef typename MappedSupernodalType::Scalar Scalar;
659  SparseLUMatrixLReturnType(const MappedSupernodalType& mapL) : m_mapL(mapL)
660  { }
661  Index rows() { return m_mapL.rows(); }
662  Index cols() { return m_mapL.cols(); }
663  template<typename Dest>
664  void solveInPlace( MatrixBase<Dest> &X) const
665  {
666    m_mapL.solveInPlace(X);
667  }
668  const MappedSupernodalType& m_mapL;
669};
670
671template<typename MatrixLType, typename MatrixUType>
672struct SparseLUMatrixUReturnType : internal::no_assignment_operator
673{
674  typedef typename MatrixLType::Index Index;
675  typedef typename MatrixLType::Scalar Scalar;
676  SparseLUMatrixUReturnType(const MatrixLType& mapL, const MatrixUType& mapU)
677  : m_mapL(mapL),m_mapU(mapU)
678  { }
679  Index rows() { return m_mapL.rows(); }
680  Index cols() { return m_mapL.cols(); }
681
682  template<typename Dest>   void solveInPlace(MatrixBase<Dest> &X) const
683  {
684    Index nrhs = X.cols();
685    Index n = X.rows();
686    // Backward solve with U
687    for (Index k = m_mapL.nsuper(); k >= 0; k--)
688    {
689      Index fsupc = m_mapL.supToCol()[k];
690      Index lda = m_mapL.colIndexPtr()[fsupc+1] - m_mapL.colIndexPtr()[fsupc]; // leading dimension
691      Index nsupc = m_mapL.supToCol()[k+1] - fsupc;
692      Index luptr = m_mapL.colIndexPtr()[fsupc];
693
694      if (nsupc == 1)
695      {
696        for (Index j = 0; j < nrhs; j++)
697        {
698          X(fsupc, j) /= m_mapL.valuePtr()[luptr];
699        }
700      }
701      else
702      {
703        Map<const Matrix<Scalar,Dynamic,Dynamic>, 0, OuterStride<> > A( &(m_mapL.valuePtr()[luptr]), nsupc, nsupc, OuterStride<>(lda) );
704        Map< Matrix<Scalar,Dynamic,Dynamic>, 0, OuterStride<> > U (&(X(fsupc,0)), nsupc, nrhs, OuterStride<>(n) );
705        U = A.template triangularView<Upper>().solve(U);
706      }
707
708      for (Index j = 0; j < nrhs; ++j)
709      {
710        for (Index jcol = fsupc; jcol < fsupc + nsupc; jcol++)
711        {
712          typename MatrixUType::InnerIterator it(m_mapU, jcol);
713          for ( ; it; ++it)
714          {
715            Index irow = it.index();
716            X(irow, j) -= X(jcol, j) * it.value();
717          }
718        }
719      }
720    } // End For U-solve
721  }
722  const MatrixLType& m_mapL;
723  const MatrixUType& m_mapU;
724};
725
726namespace internal {
727
728template<typename _MatrixType, typename Derived, typename Rhs>
729struct solve_retval<SparseLU<_MatrixType,Derived>, Rhs>
730  : solve_retval_base<SparseLU<_MatrixType,Derived>, Rhs>
731{
732  typedef SparseLU<_MatrixType,Derived> Dec;
733  EIGEN_MAKE_SOLVE_HELPERS(Dec,Rhs)
734
735  template<typename Dest> void evalTo(Dest& dst) const
736  {
737    dec()._solve(rhs(),dst);
738  }
739};
740
741template<typename _MatrixType, typename Derived, typename Rhs>
742struct sparse_solve_retval<SparseLU<_MatrixType,Derived>, Rhs>
743  : sparse_solve_retval_base<SparseLU<_MatrixType,Derived>, Rhs>
744{
745  typedef SparseLU<_MatrixType,Derived> Dec;
746  EIGEN_MAKE_SPARSE_SOLVE_HELPERS(Dec,Rhs)
747
748  template<typename Dest> void evalTo(Dest& dst) const
749  {
750    this->defaultEvalTo(dst);
751  }
752};
753} // end namespace internal
754
755} // End namespace Eigen
756
757#endif
758