SuiteSparseQRSupport.h 11.1 KB
Newer Older
Don Gagne's avatar
Don Gagne committed
1 2 3 4
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2012 Desire Nuentsa <desire.nuentsa_wakam@inria.fr>
5
// Copyright (C) 2014 Gael Guennebaud <gael.guennebaud@inria.fr>
Don Gagne's avatar
Don Gagne committed
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
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.

#ifndef EIGEN_SUITESPARSEQRSUPPORT_H
#define EIGEN_SUITESPARSEQRSUPPORT_H

namespace Eigen {
  
  template<typename MatrixType> class SPQR; 
  template<typename SPQRType> struct SPQRMatrixQReturnType; 
  template<typename SPQRType> struct SPQRMatrixQTransposeReturnType; 
  template <typename SPQRType, typename Derived> struct SPQR_QProduct;
  namespace internal {
    template <typename SPQRType> struct traits<SPQRMatrixQReturnType<SPQRType> >
    {
      typedef typename SPQRType::MatrixType ReturnType;
    };
    template <typename SPQRType> struct traits<SPQRMatrixQTransposeReturnType<SPQRType> >
    {
      typedef typename SPQRType::MatrixType ReturnType;
    };
    template <typename SPQRType, typename Derived> struct traits<SPQR_QProduct<SPQRType, Derived> >
    {
      typedef typename Derived::PlainObject ReturnType;
    };
  } // End namespace internal
  
/**
36 37 38 39 40 41 42 43 44 45 46 47 48 49 50 51 52 53 54 55 56 57 58
  * \ingroup SPQRSupport_Module
  * \class SPQR
  * \brief Sparse QR factorization based on SuiteSparseQR library
  *
  * This class is used to perform a multithreaded and multifrontal rank-revealing QR decomposition
  * of sparse matrices. The result is then used to solve linear leasts_square systems.
  * Clearly, a QR factorization is returned such that A*P = Q*R where :
  *
  * P is the column permutation. Use colsPermutation() to get it.
  *
  * Q is the orthogonal matrix represented as Householder reflectors.
  * Use matrixQ() to get an expression and matrixQ().transpose() to get the transpose.
  * You can then apply it to a vector.
  *
  * R is the sparse triangular factor. Use matrixQR() to get it as SparseMatrix.
  * NOTE : The Index type of R is always SuiteSparse_long. You can get it with SPQR::Index
  *
  * \tparam _MatrixType The type of the sparse matrix A, must be a column-major SparseMatrix<>
  *
  * \implsparsesolverconcept
  *
  *
  */
Don Gagne's avatar
Don Gagne committed
59
template<typename _MatrixType>
60
class SPQR : public SparseSolverBase<SPQR<_MatrixType> >
Don Gagne's avatar
Don Gagne committed
61
{
62 63 64
  protected:
    typedef SparseSolverBase<SPQR<_MatrixType> > Base;
    using Base::m_isInitialized;
Don Gagne's avatar
Don Gagne committed
65 66 67
  public:
    typedef typename _MatrixType::Scalar Scalar;
    typedef typename _MatrixType::RealScalar RealScalar;
68 69 70 71 72 73 74
    typedef SuiteSparse_long StorageIndex ;
    typedef SparseMatrix<Scalar, ColMajor, StorageIndex> MatrixType;
    typedef Map<PermutationMatrix<Dynamic, Dynamic, StorageIndex> > PermutationType;
    enum {
      ColsAtCompileTime = Dynamic,
      MaxColsAtCompileTime = Dynamic
    };
Don Gagne's avatar
Don Gagne committed
75 76
  public:
    SPQR() 
77
      : m_ordering(SPQR_ORDERING_DEFAULT), m_allow_tol(SPQR_DEFAULT_TOL), m_tolerance (NumTraits<Scalar>::epsilon()), m_useDefaultThreshold(true)
Don Gagne's avatar
Don Gagne committed
78 79 80 81
    { 
      cholmod_l_start(&m_cc);
    }
    
82 83
    explicit SPQR(const _MatrixType& matrix)
    : m_ordering(SPQR_ORDERING_DEFAULT), m_allow_tol(SPQR_DEFAULT_TOL), m_tolerance (NumTraits<Scalar>::epsilon()), m_useDefaultThreshold(true)
Don Gagne's avatar
Don Gagne committed
84 85 86 87 88 89 90
    {
      cholmod_l_start(&m_cc);
      compute(matrix);
    }
    
    ~SPQR()
    {
91 92 93 94 95
      SPQR_free();
      cholmod_l_finish(&m_cc);
    }
    void SPQR_free()
    {
Don Gagne's avatar
Don Gagne committed
96 97 98 99 100 101
      cholmod_l_free_sparse(&m_H, &m_cc);
      cholmod_l_free_sparse(&m_cR, &m_cc);
      cholmod_l_free_dense(&m_HTau, &m_cc);
      std::free(m_E);
      std::free(m_HPinv);
    }
102

Don Gagne's avatar
Don Gagne committed
103 104
    void compute(const _MatrixType& matrix)
    {
105 106
      if(m_isInitialized) SPQR_free();

Don Gagne's avatar
Don Gagne committed
107
      MatrixType mat(matrix);
108 109 110 111 112 113 114 115 116
      
      /* Compute the default threshold as in MatLab, see:
       * Tim Davis, "Algorithm 915, SuiteSparseQR: Multifrontal Multithreaded Rank-Revealing
       * Sparse QR Factorization, ACM Trans. on Math. Soft. 38(1), 2011, Page 8:3 
       */
      RealScalar pivotThreshold = m_tolerance;
      if(m_useDefaultThreshold) 
      {
        RealScalar max2Norm = 0.0;
117
        for (int j = 0; j < mat.cols(); j++) max2Norm = numext::maxi(max2Norm, mat.col(j).norm());
118 119 120 121
        if(max2Norm==RealScalar(0))
          max2Norm = RealScalar(1);
        pivotThreshold = 20 * (mat.rows() + mat.cols()) * max2Norm * NumTraits<RealScalar>::epsilon();
      }
Don Gagne's avatar
Don Gagne committed
122 123
      cholmod_sparse A; 
      A = viewAsCholmod(mat);
124
      m_rows = matrix.rows();
Don Gagne's avatar
Don Gagne committed
125
      Index col = matrix.cols();
126
      m_rank = SuiteSparseQR<Scalar>(m_ordering, pivotThreshold, col, &A, 
Don Gagne's avatar
Don Gagne committed
127 128 129 130
                             &m_cR, &m_E, &m_H, &m_HPinv, &m_HTau, &m_cc);

      if (!m_cR)
      {
131
        m_info = NumericalIssue;
Don Gagne's avatar
Don Gagne committed
132 133 134 135 136 137 138 139 140 141
        m_isInitialized = false;
        return;
      }
      m_info = Success;
      m_isInitialized = true;
      m_isRUpToDate = false;
    }
    /** 
     * Get the number of rows of the input matrix and the Q matrix
     */
142
    inline Index rows() const {return m_rows; }
Don Gagne's avatar
Don Gagne committed
143 144 145 146 147 148 149
    
    /** 
     * Get the number of columns of the input matrix. 
     */
    inline Index cols() const { return m_cR->ncol; }
    
    template<typename Rhs, typename Dest>
150
    void _solve_impl(const MatrixBase<Rhs> &b, MatrixBase<Dest> &dest) const
Don Gagne's avatar
Don Gagne committed
151 152 153
    {
      eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
      eigen_assert(b.cols()==1 && "This method is for vectors only");
154

Don Gagne's avatar
Don Gagne committed
155
      //Compute Q^T * b
156
      typename Dest::PlainObject y, y2;
Don Gagne's avatar
Don Gagne committed
157
      y = matrixQ().transpose() * b;
158 159
      
      // Solves with the triangular matrix R
Don Gagne's avatar
Don Gagne committed
160
      Index rk = this->rank();
161 162 163 164
      y2 = y;
      y.resize((std::max)(cols(),Index(y.rows())),y.cols());
      y.topRows(rk) = this->matrixR().topLeftCorner(rk, rk).template triangularView<Upper>().solve(y2.topRows(rk));

Don Gagne's avatar
Don Gagne committed
165
      // Apply the column permutation 
166 167 168 169 170 171 172
      // colsPermutation() performs a copy of the permutation,
      // so let's apply it manually:
      for(Index i = 0; i < rk; ++i) dest.row(m_E[i]) = y.row(i);
      for(Index i = rk; i < cols(); ++i) dest.row(m_E[i]).setZero();
      
//       y.bottomRows(y.rows()-rk).setZero();
//       dest = colsPermutation() * y.topRows(cols());
Don Gagne's avatar
Don Gagne committed
173 174 175 176 177 178 179 180 181 182
      
      m_info = Success;
    }
    
    /** \returns the sparse triangular factor R. It is a sparse matrix
     */
    const MatrixType matrixR() const
    {
      eigen_assert(m_isInitialized && " The QR factorization should be computed first, call compute()");
      if(!m_isRUpToDate) {
183
        m_R = viewAsEigen<Scalar,ColMajor, typename MatrixType::StorageIndex>(*m_cR);
Don Gagne's avatar
Don Gagne committed
184 185 186 187 188 189 190 191 192 193 194 195 196
        m_isRUpToDate = true;
      }
      return m_R;
    }
    /// Get an expression of the matrix Q
    SPQRMatrixQReturnType<SPQR> matrixQ() const
    {
      return SPQRMatrixQReturnType<SPQR>(*this);
    }
    /// Get the permutation that was applied to columns of A
    PermutationType colsPermutation() const
    { 
      eigen_assert(m_isInitialized && "Decomposition is not initialized.");
197
      return PermutationType(m_E, m_cR->ncol);
Don Gagne's avatar
Don Gagne committed
198 199 200 201 202 203 204 205 206 207 208 209 210
    }
    /**
     * Gets the rank of the matrix. 
     * It should be equal to matrixQR().cols if the matrix is full-rank
     */
    Index rank() const
    {
      eigen_assert(m_isInitialized && "Decomposition is not initialized.");
      return m_cc.SPQR_istat[4];
    }
    /// Set the fill-reducing ordering method to be used
    void setSPQROrdering(int ord) { m_ordering = ord;}
    /// Set the tolerance tol to treat columns with 2-norm < =tol as zero
211 212 213 214 215
    void setPivotThreshold(const RealScalar& tol)
    {
      m_useDefaultThreshold = false;
      m_tolerance = tol;
    }
Don Gagne's avatar
Don Gagne committed
216 217 218 219 220 221 222 223 224 225 226 227 228 229 230 231 232 233 234 235 236 237 238 239 240
    
    /** \returns a pointer to the SPQR workspace */
    cholmod_common *cholmodCommon() const { return &m_cc; }
    
    
    /** \brief Reports whether previous computation was successful.
      *
      * \returns \c Success if computation was succesful,
      *          \c NumericalIssue if the sparse QR can not be computed
      */
    ComputationInfo info() const
    {
      eigen_assert(m_isInitialized && "Decomposition is not initialized.");
      return m_info;
    }
  protected:
    bool m_analysisIsOk;
    bool m_factorizationIsOk;
    mutable bool m_isRUpToDate;
    mutable ComputationInfo m_info;
    int m_ordering; // Ordering method to use, see SPQR's manual
    int m_allow_tol; // Allow to use some tolerance during numerical factorization.
    RealScalar m_tolerance; // treat columns with 2-norm below this tolerance as zero
    mutable cholmod_sparse *m_cR; // The sparse R factor in cholmod format
    mutable MatrixType m_R; // The sparse matrix R in Eigen format
241
    mutable StorageIndex *m_E; // The permutation applied to columns
Don Gagne's avatar
Don Gagne committed
242
    mutable cholmod_sparse *m_H;  //The householder vectors
243
    mutable StorageIndex *m_HPinv; // The row permutation of H
Don Gagne's avatar
Don Gagne committed
244 245 246
    mutable cholmod_dense *m_HTau; // The Householder coefficients
    mutable Index m_rank; // The rank of the matrix
    mutable cholmod_common m_cc; // Workspace and parameters
247
    bool m_useDefaultThreshold;     // Use default threshold
248
    Index m_rows;
Don Gagne's avatar
Don Gagne committed
249 250 251 252 253 254 255
    template<typename ,typename > friend struct SPQR_QProduct;
};

template <typename SPQRType, typename Derived>
struct SPQR_QProduct : ReturnByValue<SPQR_QProduct<SPQRType,Derived> >
{
  typedef typename SPQRType::Scalar Scalar;
256
  typedef typename SPQRType::StorageIndex StorageIndex;
Don Gagne's avatar
Don Gagne committed
257 258 259 260 261 262 263 264 265 266 267 268 269 270 271 272 273 274 275 276 277 278 279 280 281 282 283 284 285 286 287 288 289 290 291 292 293 294 295 296 297 298 299 300 301 302 303 304 305 306 307 308 309 310 311 312 313
  //Define the constructor to get reference to argument types
  SPQR_QProduct(const SPQRType& spqr, const Derived& other, bool transpose) : m_spqr(spqr),m_other(other),m_transpose(transpose) {}
  
  inline Index rows() const { return m_transpose ? m_spqr.rows() : m_spqr.cols(); }
  inline Index cols() const { return m_other.cols(); }
  // Assign to a vector
  template<typename ResType>
  void evalTo(ResType& res) const
  {
    cholmod_dense y_cd;
    cholmod_dense *x_cd; 
    int method = m_transpose ? SPQR_QTX : SPQR_QX; 
    cholmod_common *cc = m_spqr.cholmodCommon();
    y_cd = viewAsCholmod(m_other.const_cast_derived());
    x_cd = SuiteSparseQR_qmult<Scalar>(method, m_spqr.m_H, m_spqr.m_HTau, m_spqr.m_HPinv, &y_cd, cc);
    res = Matrix<Scalar,ResType::RowsAtCompileTime,ResType::ColsAtCompileTime>::Map(reinterpret_cast<Scalar*>(x_cd->x), x_cd->nrow, x_cd->ncol);
    cholmod_l_free_dense(&x_cd, cc);
  }
  const SPQRType& m_spqr; 
  const Derived& m_other; 
  bool m_transpose; 
  
};
template<typename SPQRType>
struct SPQRMatrixQReturnType{
  
  SPQRMatrixQReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
  template<typename Derived>
  SPQR_QProduct<SPQRType, Derived> operator*(const MatrixBase<Derived>& other)
  {
    return SPQR_QProduct<SPQRType,Derived>(m_spqr,other.derived(),false);
  }
  SPQRMatrixQTransposeReturnType<SPQRType> adjoint() const
  {
    return SPQRMatrixQTransposeReturnType<SPQRType>(m_spqr);
  }
  // To use for operations with the transpose of Q
  SPQRMatrixQTransposeReturnType<SPQRType> transpose() const
  {
    return SPQRMatrixQTransposeReturnType<SPQRType>(m_spqr);
  }
  const SPQRType& m_spqr;
};

template<typename SPQRType>
struct SPQRMatrixQTransposeReturnType{
  SPQRMatrixQTransposeReturnType(const SPQRType& spqr) : m_spqr(spqr) {}
  template<typename Derived>
  SPQR_QProduct<SPQRType,Derived> operator*(const MatrixBase<Derived>& other)
  {
    return SPQR_QProduct<SPQRType,Derived>(m_spqr,other.derived(), true);
  }
  const SPQRType& m_spqr;
};

}// End namespace Eigen
#endif