Amesos  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
List of all members
Amesos_Klu Class Reference

Amesos_Klu: A serial, unblocked code ideal for getting started and for very sparse matrices, such as circuit matrces. More...

#include <Amesos_Umfpack.h>

Inheritance diagram for Amesos_Klu:
Inheritance graph
[legend]
Collaboration diagram for Amesos_Klu:
Collaboration graph
[legend]

Public Member Functions

 Amesos_Klu (const Epetra_LinearProblem &LinearProblem)
 Amesos_Klu Constructor. More...
 
 ~Amesos_Klu (void)
 Amesos_Klu Destructor.
 
int SymbolicFactorization ()
 Performs SymbolicFactorization on the matrix A. More...
 
int NumericFactorization ()
 Performs NumericFactorization on the matrix A. More...
 
int Solve ()
 Solves A X = B (or AT x = B) More...
 
const Epetra_LinearProblemGetProblem () const
 Get a pointer to the Problem.
 
bool MatrixShapeOK () const
 Returns true if KLU can handle this matrix shape. More...
 
int SetUseTranspose (bool UseTranspose_in)
 SetUseTranpose(true) is more efficient in Amesos_Klu. More...
 
bool UseTranspose () const
 Returns the current UseTranspose setting.
 
const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this operator.
 
int SetParameters (Teuchos::ParameterList &ParameterList)
 Updates internal variables. More...
 
int NumSymbolicFact () const
 Returns the number of symbolic factorizations performed by this object.
 
int NumNumericFact () const
 Returns the number of numeric factorizations performed by this object.
 
int NumSolve () const
 Returns the number of solves performed by this object.
 
void PrintTiming () const
 Prints timing information.
 
void PrintStatus () const
 Prints information about the factorization and solution phases.
 
void GetTiming (Teuchos::ParameterList &TimingParameterList) const
 Extracts timing information and places in parameter list.
 
- Public Member Functions inherited from Amesos_BaseSolver
virtual ~Amesos_BaseSolver ()
 Destructor.
 
virtual void setParameterList (Teuchos::RCP< Teuchos::ParameterList > const &paramList)
 Redefined from Teuchos::ParameterListAcceptor (Does Not Work)
 
virtual Teuchos::RCP
< Teuchos::ParameterList
getNonconstParameterList ()
 This is an empty stub.
 
virtual Teuchos::RCP
< Teuchos::ParameterList
unsetParameterList ()
 This is an empty stub.
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptor
virtual RCP< const ParameterListgetParameterList () const
 
virtual RCP< const ParameterListgetValidParameters () const
 

Detailed Description

Amesos_Klu: A serial, unblocked code ideal for getting started and for very sparse matrices, such as circuit matrces.

Interface to UMFPACK.

Interface to KLU internal solver.

Class Amesos_Klu is an object-oriented wrapper for KLU. KLU, whose sources are distributed within Amesos, is a serial solver for sparse matrices. KLU will solve a linear system of equations: $A X = B$, where A is an Epetra_RowMatrix and X and B are Epetra_MultiVector objects.

Amesos_Klu computes $A^T X = B$ more efficiently than $>A X = B$. The latter requires a matrix transpose – which costs both time and space.

KLU is Tim Davis' implementation of Gilbert-Peierl's left-looking sparse partial pivoting algorithm, with Eisenstat & Liu's symmetric pruning. Gilbert's version appears as [L,U,P]=lu(A) in MATLAB. It doesn't exploit dense matrix kernels, but it is the only sparse LU factorization algorithm known to be asymptotically optimal, in the sense that it takes time proportional to the number of floating-point operations. It is the precursor to SuperLU, thus the name ("clark Kent LU"). For very sparse matrices that do not suffer much fill-in (such as most circuit matrices when permuted properly) dense matrix kernels do not help, and the asymptotic run-time is of practical importance.

The klu_btf code first permutes the matrix to upper block triangular form (using two algorithms by Duff and Reid, MC13 and MC21, in the ACM Collected Algorithms). It then permutes each block via a symmetric minimum degree ordering (AMD, by Amestoy, Davis, and Duff). This ordering phase can be done just once for a sequence of matrices. Next, it factorizes each reordered block via the klu routine, which also attempts to preserve diagonal pivoting, but allows for partial pivoting if the diagonal is to small.

Date
Last updated on 24-May-05.

Constructor & Destructor Documentation

Amesos_Klu::Amesos_Klu ( const Epetra_LinearProblem LinearProblem)

Amesos_Klu Constructor.

Creates an Amesos_Klu instance, using an Epetra_LinearProblem, passing in an already-defined Epetra_LinearProblem object.

Note: The operator in LinearProblem must be an Epetra_RowMatrix.

References rcp(), and SetParameters().

Member Function Documentation

bool Amesos_Klu::MatrixShapeOK ( ) const
virtual

Returns true if KLU can handle this matrix shape.

Returns true if the matrix shape is one that KLU can handle. KLU only works with square matrices.

Implements Amesos_BaseSolver.

References GetProblem().

int Amesos_Klu::NumericFactorization ( )
virtual

Performs NumericFactorization on the matrix A.

In addition to performing numeric factorization on the matrix A, the call to NumericFactorization() implies that no change will be made to the underlying matrix without a subsequent call to NumericFactorization().

<br >Preconditions:

  • GetProblem().GetOperator() != 0 (return -1)
  • MatrixShapeOk(GetProblem().GetOperator()) == true (return -6)
  • The non-zero structure of the matrix should not have changed since the last call to SymbolicFactorization(). (return -2 if the number of non-zeros changes) Other changes can have arbitrary consequences.
  • The distribution of the matrix should not have changed since the last call to SymbolicFactorization()
  • The matrix should be indexed from 0 to n-1, unless the parameter "Reindex" was set to "true" prior to the call to SymbolicFactorization(). (return -3 - if caught)
  • The paremeter "Reindex" should not be set to "true" except on CrsMatrices. (return -4)
  • The paremeter "Reindex" should not be set to "true" unless Amesos was built with EpetraExt, i.e. with –enable-epetraext on the configure line. (return -4)
  • Internal errors retur -5.

<br >Postconditions:

  • Numeric Factorization will be performed (or marked to be performed) allowing Solve() to be performed correctly despite a potential change in in the matrix values (though not in the non-zero structure).
Returns
Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

References Amesos_Time::AddTime(), Amesos_Control::AddZeroToDiag_, Epetra_LinearProblem::GetMatrix(), Amesos_Status::IsNumericFactorizationOK_, Amesos_Status::IsSymbolicFactorizationOK_, Amesos_Status::NumNumericFact_, Amesos_Time::ResetTimer(), and SymbolicFactorization().

Referenced by Solve().

int Amesos_Klu::SetParameters ( Teuchos::ParameterList ParameterList)
virtual

Updates internal variables.

  <br \>Preconditions:<ul>
  <li>None.</li>
  </ul>

  <br \>Postconditions:<ul> 
  <li>Internal variables controlling the factorization and solve will
  be updated and take effect on all subseuent calls to NumericFactorization() 
  and Solve().</li>
  <li>All parameters whose value are to differ from the default values must 

be included in ParameterList. Parameters not specified in ParameterList revert to their default values.

Returns
Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

References Teuchos::ParameterList::get(), Teuchos::ParameterList::isParameter(), Teuchos::ParameterList::isSublist(), and Teuchos::ParameterList::sublist().

Referenced by Amesos_Klu().

int Amesos_Klu::SetUseTranspose ( bool  UseTranspose_in)
inlinevirtual

SetUseTranpose(true) is more efficient in Amesos_Klu.

If SetUseTranspose() is set to true, $A^T X = B$ is computed.

Implements Amesos_BaseSolver.

int Amesos_Klu::Solve ( )
virtual

Solves A X = B (or AT x = B)

<br >Preconditions:

<br >Postconditions:

  • X will be set such that A X = B (or AT X = B), within the limits of the accuracy of the underlying solver.
Returns
Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

References Amesos_Time::AddTime(), Amesos_Utils::ComputeTrueResidual(), Amesos_Status::ComputeTrueResidual_, Amesos_Utils::ComputeVectorNorms(), Amesos_Status::ComputeVectorNorms_, Epetra_LinearProblem::GetLHS(), Epetra_LinearProblem::GetMatrix(), Epetra_LinearProblem::GetRHS(), Insert, Amesos_Status::IsNumericFactorizationOK_, NumericFactorization(), Amesos_Status::NumSolve_, rcp(), Teuchos::rcp(), Amesos_Control::Reindex_, Amesos_Time::ResetTimer(), and UseTranspose().

int Amesos_Klu::SymbolicFactorization ( )
virtual

Performs SymbolicFactorization on the matrix A.

In addition to performing symbolic factorization on the matrix A, the call to SymbolicFactorization() implies that no change will be made to the non-zero structure of the underlying matrix without a subsequent call to SymbolicFactorization().

<br >Preconditions:

<br >Postconditions:

Returns
Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

References Amesos_Time::AddTime(), Comm(), Amesos_Time::CreateTimer(), Epetra_LinearProblem::GetLHS(), Epetra_LinearProblem::GetMatrix(), Epetra_LinearProblem::GetRHS(), Amesos_Status::IsNumericFactorizationOK_, Amesos_Status::IsSymbolicFactorizationOK_, Epetra_Comm::MyPID(), Epetra_Comm::NumProc(), Amesos_Status::NumSymbolicFact_, Teuchos::rcp(), Amesos_Control::Reindex_, and Amesos_Time::ResetTimer().

Referenced by NumericFactorization().


The documentation for this class was generated from the following files: