MOAB  4.9.3pre
Eigen::UmfPackLU< _MatrixType > Class Template Reference

A sparse LU factorization and solver based on UmfPack. More...

#include <UmfPackSupport.h>

Inheritance diagram for Eigen::UmfPackLU< _MatrixType >:
Collaboration diagram for Eigen::UmfPackLU< _MatrixType >:

List of all members.

Public Types

enum  { ColsAtCompileTime = MatrixType::ColsAtCompileTime, MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime }
typedef _MatrixType MatrixType
typedef MatrixType::Scalar Scalar
typedef MatrixType::RealScalar RealScalar
typedef MatrixType::StorageIndex StorageIndex
typedef Matrix< Scalar,
Dynamic, 1 > 
Vector
typedef Matrix< int,
1, MatrixType::ColsAtCompileTime > 
IntRowVectorType
typedef Matrix< int,
MatrixType::RowsAtCompileTime, 1 > 
IntColVectorType
typedef SparseMatrix< ScalarLUMatrixType
typedef SparseMatrix< Scalar,
ColMajor, int > 
UmfpackMatrixType
typedef Ref< const
UmfpackMatrixType,
StandardCompressedFormat
UmfpackMatrixRef
typedef Array< double,
UMFPACK_CONTROL, 1 > 
UmfpackControl

Public Member Functions

 UmfPackLU ()
template<typename InputMatrixType >
 UmfPackLU (const InputMatrixType &matrix)
 ~UmfPackLU ()
Index rows () const
Index cols () const
ComputationInfo info () const
 Reports whether previous computation was successful.
const LUMatrixTypematrixL () const
const LUMatrixTypematrixU () const
const IntColVectorTypepermutationP () const
const IntRowVectorTypepermutationQ () const
template<typename InputMatrixType >
void compute (const InputMatrixType &matrix)
template<typename InputMatrixType >
void analyzePattern (const InputMatrixType &matrix)
int umfpackFactorizeReturncode () const
const UmfpackControlumfpackControl () const
UmfpackControlumfpackControl ()
template<typename InputMatrixType >
void factorize (const InputMatrixType &matrix)
template<typename BDerived , typename XDerived >
bool _solve_impl (const MatrixBase< BDerived > &b, MatrixBase< XDerived > &x) const
Scalar determinant () const
void extractData () const

Protected Types

typedef SparseSolverBase
< UmfPackLU< _MatrixType > > 
Base

Protected Member Functions

void init ()
void analyzePattern_impl ()
void factorize_impl ()
template<typename MatrixDerived >
void grab (const EigenBase< MatrixDerived > &A)
void grab (const UmfpackMatrixRef &A)

Protected Attributes

LUMatrixType m_l
int m_fact_errorCode
UmfpackControl m_control
LUMatrixType m_u
IntColVectorType m_p
IntRowVectorType m_q
UmfpackMatrixType m_dummy
UmfpackMatrixRef mp_matrix
void * m_numeric
void * m_symbolic
ComputationInfo m_info
int m_factorizationIsOk
int m_analysisIsOk
bool m_extractedDataAreDirty

Private Member Functions

 UmfPackLU (UmfPackLU &)

Detailed Description

template<typename _MatrixType>
class Eigen::UmfPackLU< _MatrixType >

A sparse LU factorization and solver based on UmfPack.

This class allows to solve for A.X = B sparse linear problems via a LU factorization using the UmfPack library. The sparse matrix A must be squared and full rank. The vectors or matrices X and B can be either dense or sparse.

Warning:
The input matrix A should be in a compressed and column-major form. Otherwise an expensive copy will be made. You can call the inexpensive makeCompressed() to get a compressed matrix.
Template Parameters:
_MatrixTypethe type of the sparse matrix A, it must be a SparseMatrix<>
See also:
TutorialSparseSolverConcept, class SparseLU

Definition at line 134 of file UmfPackSupport.h.


Member Typedef Documentation

template<typename _MatrixType>
typedef SparseSolverBase<UmfPackLU<_MatrixType> > Eigen::UmfPackLU< _MatrixType >::Base [protected]

Definition at line 137 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef Matrix<int, MatrixType::RowsAtCompileTime, 1> Eigen::UmfPackLU< _MatrixType >::IntColVectorType

Definition at line 147 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef Matrix<int, 1, MatrixType::ColsAtCompileTime> Eigen::UmfPackLU< _MatrixType >::IntRowVectorType

Definition at line 146 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef SparseMatrix<Scalar> Eigen::UmfPackLU< _MatrixType >::LUMatrixType

Definition at line 148 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef _MatrixType Eigen::UmfPackLU< _MatrixType >::MatrixType

Definition at line 141 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef MatrixType::RealScalar Eigen::UmfPackLU< _MatrixType >::RealScalar

Definition at line 143 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef MatrixType::Scalar Eigen::UmfPackLU< _MatrixType >::Scalar

Definition at line 142 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef MatrixType::StorageIndex Eigen::UmfPackLU< _MatrixType >::StorageIndex

Definition at line 144 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef Array<double, UMFPACK_CONTROL, 1> Eigen::UmfPackLU< _MatrixType >::UmfpackControl

Definition at line 158 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef Ref<const UmfpackMatrixType, StandardCompressedFormat> Eigen::UmfPackLU< _MatrixType >::UmfpackMatrixRef

Definition at line 150 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef SparseMatrix<Scalar,ColMajor,int> Eigen::UmfPackLU< _MatrixType >::UmfpackMatrixType

Definition at line 149 of file UmfPackSupport.h.

template<typename _MatrixType>
typedef Matrix<Scalar,Dynamic,1> Eigen::UmfPackLU< _MatrixType >::Vector

Definition at line 145 of file UmfPackSupport.h.


Member Enumeration Documentation

template<typename _MatrixType>
anonymous enum
Enumerator:
ColsAtCompileTime 
MaxColsAtCompileTime 

Definition at line 151 of file UmfPackSupport.h.

         {
      ColsAtCompileTime = MatrixType::ColsAtCompileTime,
      MaxColsAtCompileTime = MatrixType::MaxColsAtCompileTime
    };

Constructor & Destructor Documentation

template<typename _MatrixType>
Eigen::UmfPackLU< _MatrixType >::UmfPackLU ( ) [inline]

Definition at line 160 of file UmfPackSupport.h.

      : m_dummy(0,0), mp_matrix(m_dummy)
    {
      init();
    }
template<typename _MatrixType>
template<typename InputMatrixType >
Eigen::UmfPackLU< _MatrixType >::UmfPackLU ( const InputMatrixType &  matrix) [inline, explicit]

Definition at line 167 of file UmfPackSupport.h.

      : mp_matrix(matrix)
    {
      init();
      compute(matrix);
    }
template<typename _MatrixType>
Eigen::UmfPackLU< _MatrixType >::~UmfPackLU ( ) [inline]
template<typename _MatrixType>
Eigen::UmfPackLU< _MatrixType >::UmfPackLU ( UmfPackLU< _MatrixType > &  ) [inline, private]

Definition at line 382 of file UmfPackSupport.h.

{ }

Member Function Documentation

template<typename MatrixType >
template<typename BDerived , typename XDerived >
bool Eigen::UmfPackLU< MatrixType >::_solve_impl ( const MatrixBase< BDerived > &  b,
MatrixBase< XDerived > &  x 
) const

Definition at line 424 of file UmfPackSupport.h.

{
  Index rhsCols = b.cols();
  eigen_assert((BDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major rhs yet");
  eigen_assert((XDerived::Flags&RowMajorBit)==0 && "UmfPackLU backend does not support non col-major result yet");
  eigen_assert(b.derived().data() != x.derived().data() && " Umfpack does not support inplace solve");
  
  int errorCode;
  Scalar* x_ptr = 0;
  Matrix<Scalar,Dynamic,1> x_tmp;
  if(x.innerStride()!=1)
  {
    x_tmp.resize(x.rows());
    x_ptr = x_tmp.data();
  }
  for (int j=0; j<rhsCols; ++j)
  {
    if(x.innerStride()==1)
      x_ptr = &x.col(j).coeffRef(0);
    errorCode = umfpack_solve(UMFPACK_A,
        mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
        x_ptr, &b.const_cast_derived().col(j).coeffRef(0), m_numeric, m_control.data(), 0);
    if(x.innerStride()!=1)
      x.col(j) = x_tmp;
    if (errorCode!=0)
      return false;
  }

  return true;
}
template<typename _MatrixType>
template<typename InputMatrixType >
void Eigen::UmfPackLU< _MatrixType >::analyzePattern ( const InputMatrixType &  matrix) [inline]

Performs a symbolic decomposition on the sparcity of matrix.

This function is particularly useful when solving for several problems having the same structure.

See also:
factorize(), compute()

Definition at line 239 of file UmfPackSupport.h.

template<typename _MatrixType>
void Eigen::UmfPackLU< _MatrixType >::analyzePattern_impl ( ) [inline, protected]

Definition at line 319 of file UmfPackSupport.h.

    {
      umfpack_defaults(m_control.data(), Scalar());
      int errorCode = 0;
      errorCode = umfpack_symbolic(internal::convert_index<int>(mp_matrix.rows()),
                                   internal::convert_index<int>(mp_matrix.cols()),
                                   mp_matrix.outerIndexPtr(), mp_matrix.innerIndexPtr(), mp_matrix.valuePtr(),
                                   &m_symbolic, m_control.data(), 0);

      m_isInitialized = true;
      m_info = errorCode ? InvalidInput : Success;
      m_analysisIsOk = true;
      m_factorizationIsOk = false;
      m_extractedDataAreDirty = true;
    }
template<typename _MatrixType>
Index Eigen::UmfPackLU< _MatrixType >::cols ( void  ) const [inline]

Definition at line 181 of file UmfPackSupport.h.

{ return mp_matrix.cols(); }
template<typename _MatrixType>
template<typename InputMatrixType >
void Eigen::UmfPackLU< _MatrixType >::compute ( const InputMatrixType &  matrix) [inline]

Computes the sparse Cholesky decomposition of matrix Note that the matrix should be column-major, and in compressed format for best performance.

See also:
SparseMatrix::makeCompressed().

Definition at line 223 of file UmfPackSupport.h.

template<typename MatrixType >
UmfPackLU< MatrixType >::Scalar Eigen::UmfPackLU< MatrixType >::determinant ( ) const

Definition at line 415 of file UmfPackSupport.h.

{
  Scalar det;
  umfpack_get_determinant(&det, 0, m_numeric, 0);
  return det;
}
template<typename MatrixType >
void Eigen::UmfPackLU< MatrixType >::extractData ( ) const

Definition at line 387 of file UmfPackSupport.h.

{
  if (m_extractedDataAreDirty)
  {
    // get size of the data
    int lnz, unz, rows, cols, nz_udiag;
    umfpack_get_lunz(&lnz, &unz, &rows, &cols, &nz_udiag, m_numeric, Scalar());

    // allocate data
    m_l.resize(rows,(std::min)(rows,cols));
    m_l.resizeNonZeros(lnz);

    m_u.resize((std::min)(rows,cols),cols);
    m_u.resizeNonZeros(unz);

    m_p.resize(rows);
    m_q.resize(cols);

    // extract
    umfpack_get_numeric(m_l.outerIndexPtr(), m_l.innerIndexPtr(), m_l.valuePtr(),
                        m_u.outerIndexPtr(), m_u.innerIndexPtr(), m_u.valuePtr(),
                        m_p.data(), m_q.data(), 0, 0, 0, m_numeric);

    m_extractedDataAreDirty = false;
  }
}
template<typename _MatrixType>
template<typename InputMatrixType >
void Eigen::UmfPackLU< _MatrixType >::factorize ( const InputMatrixType &  matrix) [inline]

Performs a numeric decomposition of matrix

The given matrix must has the same sparcity than the matrix on which the pattern anylysis has been performed.

See also:
analyzePattern(), compute()

Definition at line 289 of file UmfPackSupport.h.

    {
      eigen_assert(m_analysisIsOk && "UmfPackLU: you must first call analyzePattern()");
      if(m_numeric)
        umfpack_free_numeric(&m_numeric,Scalar());

      grab(matrix.derived());
      
      factorize_impl();
    }
template<typename _MatrixType>
void Eigen::UmfPackLU< _MatrixType >::factorize_impl ( ) [inline, protected]

Definition at line 335 of file UmfPackSupport.h.

template<typename _MatrixType>
template<typename MatrixDerived >
void Eigen::UmfPackLU< _MatrixType >::grab ( const EigenBase< MatrixDerived > &  A) [inline, protected]

Definition at line 346 of file UmfPackSupport.h.

    {
      mp_matrix.~UmfpackMatrixRef();
      ::new (&mp_matrix) UmfpackMatrixRef(A.derived());
    }
template<typename _MatrixType>
void Eigen::UmfPackLU< _MatrixType >::grab ( const UmfpackMatrixRef A) [inline, protected]

Definition at line 352 of file UmfPackSupport.h.

    {
      if(&(A.derived()) != &mp_matrix)
      {
        mp_matrix.~UmfpackMatrixRef();
        ::new (&mp_matrix) UmfpackMatrixRef(A);
      }
    }
template<typename _MatrixType>
ComputationInfo Eigen::UmfPackLU< _MatrixType >::info ( ) const [inline]

Reports whether previous computation was successful.

Returns:
Success if computation was succesful, NumericalIssue if the matrix.appears to be negative.

Definition at line 188 of file UmfPackSupport.h.

    {
      eigen_assert(m_isInitialized && "Decomposition is not initialized.");
      return m_info;
    }
template<typename _MatrixType>
void Eigen::UmfPackLU< _MatrixType >::init ( ) [inline, protected]

Definition at line 310 of file UmfPackSupport.h.

template<typename _MatrixType>
const LUMatrixType& Eigen::UmfPackLU< _MatrixType >::matrixL ( ) const [inline]

Definition at line 194 of file UmfPackSupport.h.

template<typename _MatrixType>
const LUMatrixType& Eigen::UmfPackLU< _MatrixType >::matrixU ( ) const [inline]

Definition at line 200 of file UmfPackSupport.h.

template<typename _MatrixType>
const IntColVectorType& Eigen::UmfPackLU< _MatrixType >::permutationP ( ) const [inline]

Definition at line 206 of file UmfPackSupport.h.

template<typename _MatrixType>
const IntRowVectorType& Eigen::UmfPackLU< _MatrixType >::permutationQ ( ) const [inline]

Definition at line 212 of file UmfPackSupport.h.

template<typename _MatrixType>
Index Eigen::UmfPackLU< _MatrixType >::rows ( void  ) const [inline]

Definition at line 180 of file UmfPackSupport.h.

{ return mp_matrix.rows(); }
template<typename _MatrixType>
const UmfpackControl& Eigen::UmfPackLU< _MatrixType >::umfpackControl ( ) const [inline]

Provides access to the control settings array used by UmfPack.

If this array contains NaN's, the default values are used.

See UMFPACK documentation for details.

Definition at line 266 of file UmfPackSupport.h.

    {
      return m_control;
    }
template<typename _MatrixType>
UmfpackControl& Eigen::UmfPackLU< _MatrixType >::umfpackControl ( ) [inline]

Provides access to the control settings array used by UmfPack.

If this array contains NaN's, the default values are used.

See UMFPACK documentation for details.

Definition at line 277 of file UmfPackSupport.h.

    {
      return m_control;
    }
template<typename _MatrixType>
int Eigen::UmfPackLU< _MatrixType >::umfpackFactorizeReturncode ( ) const [inline]

Provides the return status code returned by UmfPack during the numeric factorization.

See also:
factorize(), compute()

Definition at line 254 of file UmfPackSupport.h.

    {
      eigen_assert(m_numeric && "UmfPackLU: you must first call factorize()");
      return m_fact_errorCode;
    }

Member Data Documentation

template<typename _MatrixType>
int Eigen::UmfPackLU< _MatrixType >::m_analysisIsOk [protected]

Definition at line 378 of file UmfPackSupport.h.

template<typename _MatrixType>
UmfpackControl Eigen::UmfPackLU< _MatrixType >::m_control [protected]

Definition at line 364 of file UmfPackSupport.h.

template<typename _MatrixType>
UmfpackMatrixType Eigen::UmfPackLU< _MatrixType >::m_dummy [protected]

Definition at line 370 of file UmfPackSupport.h.

template<typename _MatrixType>
bool Eigen::UmfPackLU< _MatrixType >::m_extractedDataAreDirty [mutable, protected]

Definition at line 379 of file UmfPackSupport.h.

template<typename _MatrixType>
int Eigen::UmfPackLU< _MatrixType >::m_fact_errorCode [protected]

Definition at line 363 of file UmfPackSupport.h.

template<typename _MatrixType>
int Eigen::UmfPackLU< _MatrixType >::m_factorizationIsOk [protected]

Definition at line 377 of file UmfPackSupport.h.

template<typename _MatrixType>
ComputationInfo Eigen::UmfPackLU< _MatrixType >::m_info [mutable, protected]

Definition at line 376 of file UmfPackSupport.h.

template<typename _MatrixType>
LUMatrixType Eigen::UmfPackLU< _MatrixType >::m_l [mutable, protected]

Definition at line 362 of file UmfPackSupport.h.

template<typename _MatrixType>
void* Eigen::UmfPackLU< _MatrixType >::m_numeric [protected]

Definition at line 373 of file UmfPackSupport.h.

template<typename _MatrixType>
IntColVectorType Eigen::UmfPackLU< _MatrixType >::m_p [mutable, protected]

Definition at line 367 of file UmfPackSupport.h.

template<typename _MatrixType>
IntRowVectorType Eigen::UmfPackLU< _MatrixType >::m_q [mutable, protected]

Definition at line 368 of file UmfPackSupport.h.

template<typename _MatrixType>
void* Eigen::UmfPackLU< _MatrixType >::m_symbolic [protected]

Definition at line 374 of file UmfPackSupport.h.

template<typename _MatrixType>
LUMatrixType Eigen::UmfPackLU< _MatrixType >::m_u [mutable, protected]

Definition at line 366 of file UmfPackSupport.h.

template<typename _MatrixType>
UmfpackMatrixRef Eigen::UmfPackLU< _MatrixType >::mp_matrix [protected]

Definition at line 371 of file UmfPackSupport.h.


The documentation for this class was generated from the following file:
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Defines