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

A sparse direct LU factorization and solver based on the SuperLU library. More...

#include <SuperLUSupport.h>

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

List of all members.

Public Types

typedef SuperLUBase
< _MatrixType, SuperLU
Base
typedef _MatrixType MatrixType
typedef Base::Scalar Scalar
typedef Base::RealScalar RealScalar
typedef Base::StorageIndex StorageIndex
typedef Base::IntRowVectorType IntRowVectorType
typedef Base::IntColVectorType IntColVectorType
typedef Base::PermutationMap PermutationMap
typedef Base::LUMatrixType LUMatrixType
typedef TriangularView
< LUMatrixType, Lower|UnitDiag
LMatrixType
typedef TriangularView
< LUMatrixType, Upper
UMatrixType

Public Member Functions

 SuperLU ()
 SuperLU (const MatrixType &matrix)
 ~SuperLU ()
void analyzePattern (const MatrixType &matrix)
void factorize (const MatrixType &matrix)
template<typename Rhs , typename Dest >
void _solve_impl (const MatrixBase< Rhs > &b, MatrixBase< Dest > &dest) const
const LMatrixTypematrixL () const
const UMatrixTypematrixU () const
const IntColVectorTypepermutationP () const
const IntRowVectorTypepermutationQ () const
Scalar determinant () const

Protected Member Functions

void init ()

Private Member Functions

 SuperLU (SuperLU &)

Detailed Description

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

A sparse direct LU factorization and solver based on the SuperLU library.

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

Template Parameters:
_MatrixTypethe type of the sparse matrix A, it must be a SparseMatrix<>
Warning:
This class is only for the 4.x versions of SuperLU. The 3.x and 5.x versions are not supported.
See also:
TutorialSparseSolverConcept, class SparseLU

Definition at line 462 of file SuperLUSupport.h.


Member Typedef Documentation

template<typename _MatrixType>
typedef SuperLUBase<_MatrixType,SuperLU> Eigen::SuperLU< _MatrixType >::Base

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 465 of file SuperLUSupport.h.

template<typename _MatrixType>
typedef Base::IntColVectorType Eigen::SuperLU< _MatrixType >::IntColVectorType

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 471 of file SuperLUSupport.h.

template<typename _MatrixType>
typedef Base::IntRowVectorType Eigen::SuperLU< _MatrixType >::IntRowVectorType

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 470 of file SuperLUSupport.h.

template<typename _MatrixType>
typedef TriangularView<LUMatrixType, Lower|UnitDiag> Eigen::SuperLU< _MatrixType >::LMatrixType

Definition at line 474 of file SuperLUSupport.h.

template<typename _MatrixType>
typedef Base::LUMatrixType Eigen::SuperLU< _MatrixType >::LUMatrixType

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 473 of file SuperLUSupport.h.

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

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 466 of file SuperLUSupport.h.

template<typename _MatrixType>
typedef Base::PermutationMap Eigen::SuperLU< _MatrixType >::PermutationMap

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 472 of file SuperLUSupport.h.

template<typename _MatrixType>
typedef Base::RealScalar Eigen::SuperLU< _MatrixType >::RealScalar

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 468 of file SuperLUSupport.h.

template<typename _MatrixType>
typedef Base::Scalar Eigen::SuperLU< _MatrixType >::Scalar

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 467 of file SuperLUSupport.h.

template<typename _MatrixType>
typedef Base::StorageIndex Eigen::SuperLU< _MatrixType >::StorageIndex

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 469 of file SuperLUSupport.h.

template<typename _MatrixType>
typedef TriangularView<LUMatrixType, Upper> Eigen::SuperLU< _MatrixType >::UMatrixType

Definition at line 475 of file SuperLUSupport.h.


Constructor & Destructor Documentation

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

Definition at line 480 of file SuperLUSupport.h.

: Base() { init(); }
template<typename _MatrixType>
Eigen::SuperLU< _MatrixType >::SuperLU ( const MatrixType matrix) [inline, explicit]

Definition at line 482 of file SuperLUSupport.h.

                                               : Base()
    {
      init();
      Base::compute(matrix);
    }
template<typename _MatrixType>
Eigen::SuperLU< _MatrixType >::~SuperLU ( ) [inline]

Definition at line 488 of file SuperLUSupport.h.

    {
    }
template<typename _MatrixType>
Eigen::SuperLU< _MatrixType >::SuperLU ( SuperLU< _MatrixType > &  ) [inline, private]

Definition at line 583 of file SuperLUSupport.h.

{ }

Member Function Documentation

template<typename MatrixType >
template<typename Rhs , typename Dest >
void Eigen::SuperLU< MatrixType >::_solve_impl ( const MatrixBase< Rhs > &  b,
MatrixBase< Dest > &  dest 
) const

Definition at line 623 of file SuperLUSupport.h.

{
  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for solving, you must first call either compute() or analyzePattern()/factorize()");

  const Index size = m_matrix.rows();
  const Index rhsCols = b.cols();
  eigen_assert(size==b.rows());

  m_sluOptions.Trans = NOTRANS;
  m_sluOptions.Fact = FACTORED;
  m_sluOptions.IterRefine = NOREFINE;
  

  m_sluFerr.resize(rhsCols);
  m_sluBerr.resize(rhsCols);
  
  Ref<const Matrix<typename Rhs::Scalar,Dynamic,Dynamic,ColMajor> > b_ref(b);
  Ref<const Matrix<typename Dest::Scalar,Dynamic,Dynamic,ColMajor> > x_ref(x);
  
  m_sluB = SluMatrix::Map(b_ref.const_cast_derived());
  m_sluX = SluMatrix::Map(x_ref.const_cast_derived());
  
  typename Rhs::PlainObject b_cpy;
  if(m_sluEqued!='N')
  {
    b_cpy = b;
    m_sluB = SluMatrix::Map(b_cpy.const_cast_derived());  
  }

  StatInit(&m_sluStat);
  int info = 0;
  RealScalar recip_pivot_growth, rcond;
  SuperLU_gssvx(&m_sluOptions, &m_sluA,
                m_q.data(), m_p.data(),
                &m_sluEtree[0], &m_sluEqued,
                &m_sluRscale[0], &m_sluCscale[0],
                &m_sluL, &m_sluU,
                NULL, 0,
                &m_sluB, &m_sluX,
                &recip_pivot_growth, &rcond,
                &m_sluFerr[0], &m_sluBerr[0],
                &m_sluStat, &info, Scalar());
  StatFree(&m_sluStat);
  
  if(x.derived().data() != x_ref.data())
    x = x_ref;
  
  m_info = info==0 ? Success : NumericalIssue;
}
template<typename _MatrixType>
void Eigen::SuperLU< _MatrixType >::analyzePattern ( const MatrixType 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()

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 498 of file SuperLUSupport.h.

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

Definition at line 767 of file SuperLUSupport.h.

{
  eigen_assert(m_factorizationIsOk && "The decomposition is not in a valid state for computing the determinant, you must first call either compute() or analyzePattern()/factorize()");
  
  if (m_extractedDataAreDirty)
    this->extractData();

  Scalar det = Scalar(1);
  for (int j=0; j<m_u.cols(); ++j)
  {
    if (m_u.outerIndexPtr()[j+1]-m_u.outerIndexPtr()[j] > 0)
    {
      int lastId = m_u.outerIndexPtr()[j+1]-1;
      eigen_assert(m_u.innerIndexPtr()[lastId]<=j);
      if (m_u.innerIndexPtr()[lastId]==j)
        det *= m_u.valuePtr()[lastId];
    }
  }
  if(PermutationMap(m_p.data(),m_p.size()).determinant()*PermutationMap(m_q.data(),m_q.size()).determinant()<0)
    det = -det;
  if(m_sluEqued!='N')
    return det/m_sluRscale.prod()/m_sluCscale.prod();
  else
    return det;
}
template<typename MatrixType >
void Eigen::SuperLU< MatrixType >::factorize ( const MatrixType matrix)

Performs a numeric decomposition of matrix

The given matrix must has the same sparcity than the matrix on which the symbolic decomposition has been performed.

See also:
analyzePattern()

Definition at line 587 of file SuperLUSupport.h.

{
  eigen_assert(m_analysisIsOk && "You must first call analyzePattern()");
  if(!m_analysisIsOk)
  {
    m_info = InvalidInput;
    return;
  }
  
  this->initFactorization(a);
  
  m_sluOptions.ColPerm = COLAMD;
  int info = 0;
  RealScalar recip_pivot_growth, rcond;
  RealScalar ferr, berr;

  StatInit(&m_sluStat);
  SuperLU_gssvx(&m_sluOptions, &m_sluA, m_q.data(), m_p.data(), &m_sluEtree[0],
                &m_sluEqued, &m_sluRscale[0], &m_sluCscale[0],
                &m_sluL, &m_sluU,
                NULL, 0,
                &m_sluB, &m_sluX,
                &recip_pivot_growth, &rcond,
                &ferr, &berr,
                &m_sluStat, &info, Scalar());
  StatFree(&m_sluStat);

  m_extractedDataAreDirty = true;

  // FIXME how to better check for errors ???
  m_info = info == 0 ? Success : NumericalIssue;
  m_factorizationIsOk = true;
}
template<typename _MatrixType>
void Eigen::SuperLU< _MatrixType >::init ( ) [inline, protected]

Reimplemented from Eigen::SuperLUBase< _MatrixType, SuperLU< _MatrixType > >.

Definition at line 570 of file SuperLUSupport.h.

    {
      Base::init();
      
      set_default_options(&this->m_sluOptions);
      m_sluOptions.PrintStat        = NO;
      m_sluOptions.ConditionNumber  = NO;
      m_sluOptions.Trans            = NOTRANS;
      m_sluOptions.ColPerm          = COLAMD;
    }
template<typename _MatrixType>
const LMatrixType& Eigen::SuperLU< _MatrixType >::matrixL ( ) const [inline]

Definition at line 517 of file SuperLUSupport.h.

    {
      if (m_extractedDataAreDirty) this->extractData();
      return m_l;
    }
template<typename _MatrixType>
const UMatrixType& Eigen::SuperLU< _MatrixType >::matrixU ( ) const [inline]

Definition at line 523 of file SuperLUSupport.h.

    {
      if (m_extractedDataAreDirty) this->extractData();
      return m_u;
    }
template<typename _MatrixType>
const IntColVectorType& Eigen::SuperLU< _MatrixType >::permutationP ( ) const [inline]

Definition at line 529 of file SuperLUSupport.h.

    {
      if (m_extractedDataAreDirty) this->extractData();
      return m_p;
    }
template<typename _MatrixType>
const IntRowVectorType& Eigen::SuperLU< _MatrixType >::permutationQ ( ) const [inline]

Definition at line 535 of file SuperLUSupport.h.

    {
      if (m_extractedDataAreDirty) this->extractData();
      return m_q;
    }

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