Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
Amesos_BaseSolver Class Referenceabstract

Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators. More...

#include <Amesos_BaseSolver.h>

Inheritance diagram for Amesos_BaseSolver:
Inheritance graph
[legend]
virtual ~Amesos_BaseSolver ()
 Destructor. More...
 
virtual int SymbolicFactorization ()=0
 Performs SymbolicFactorization on the matrix A. More...
 
virtual int NumericFactorization ()=0
 Performs NumericFactorization on the matrix A. More...
 
virtual int Solve ()=0
 Solves A X = B (or AT x = B) More...
 
virtual int SetUseTranspose (bool UseTranspose)=0
 If set true, X will be set to the solution of AT X = B (not A X = B) More...
 
virtual bool UseTranspose () const =0
 Returns the current UseTranspose setting. More...
 
virtual int SetParameters (Teuchos::ParameterList &ParameterList)=0
 Updates internal variables. More...
 
virtual const
Epetra_LinearProblem
GetProblem () const =0
 Returns the Epetra_LinearProblem. More...
 
virtual bool MatrixShapeOK () const =0
 Returns true if the solver can handle this matrix shape. More...
 
virtual const Epetra_CommComm () const =0
 Returns a pointer to the Epetra_Comm communicator associated with this operator. More...
 
virtual int NumSymbolicFact () const =0
 Returns the number of symbolic factorizations performed by this object. More...
 
virtual int NumNumericFact () const =0
 Returns the number of numeric factorizations performed by this object. More...
 
virtual int NumSolve () const =0
 Returns the number of solves performed by this object. More...
 
virtual void PrintStatus () const =0
 Prints status information about the current solver. More...
 
virtual void PrintTiming () const =0
 Prints timing information about the current solver. More...
 
virtual void setParameterList (Teuchos::RCP< Teuchos::ParameterList > const &paramList)
 Redefined from Teuchos::ParameterListAcceptor (Does Not Work) More...
 
virtual Teuchos::RCP
< Teuchos::ParameterList
getNonconstParameterList ()
 This is an empty stub. More...
 
virtual Teuchos::RCP
< Teuchos::ParameterList
unsetParameterList ()
 This is an empty stub. More...
 
virtual void GetTiming (Teuchos::ParameterList &TimingParameterList) const
 Extracts timing information from the current solver and places it in the parameter list. (Does Not Work) More...
 

Additional Inherited Members

- Public Member Functions inherited from Teuchos::ParameterListAcceptor
virtual RCP< const ParameterListgetParameterList () const
 
virtual RCP< const ParameterListgetValidParameters () const
 

Detailed Description

Amesos_BaseSolver: A pure virtual class for direct solution of real-valued double-precision operators.

Pure virtual class for all Amesos concrete implementions.

The Amesos_BaseSolver class is a pure virtual class (that is, it specifies interface only) that enables the use of real-valued double-precision direct sparse solvers. Every Amesos class named Amesos_SolverName derives from Amesos_BaseSolver.

Usage Examples

Basic calling sequence

The basic calling sequence solves A x = b or AT x = b without specifying how A has changed between each call to Solve().

Amesos_SolverName Solver(Problem);
Problem.SetUseTranspose(false);
while( ... ) {
<Here code which may change A or B>
Solver.SymbolicFactorization();
Solver.NumericFactorization();
Solver.Solve();
}

Re-using the symbolic factorization

The following calling sequence performs multiple solves of A x = b or AT x = b in cases where the non-zero structure of A remains unchanged between each call to Solve().

Amesos_SolverName Solver(Problem);
Problem.SetUseTranspose(false);
Solver.SymbolicFactorization();
while( ... ) {
<Here code which may change B or the non-zero values
(but not the non-zero structure) of A>
Solver.NumericFactorization() ;
Solver.Solve() ;
}

Re-using the numeric factorization

The following calling sequence performs multiple solves of A x = b or AT x = b provided that A remains unchanged between each call to Solve().

Amesos_SolverName Solver(Problem);
Problem.SetUseTranspose( false );
Solver.NumericFactorization();
while( ... ) {
<Here code which may change B but not A>
Solver.Solve();
}

Constructor requirements

Every Amesos_SolverName class should accept an Epetra_LinearProblem

Mathematical methods

Four mathematical methods are defined in the base class Amesos_BaseSolver: SymbolicFactorization(), NumericFactorization(), and Solve().

Switching concrete classes

Different concrete classes, each based on a different third party solver, will have different performance characteristics and will accept different parameters.

Changing the values of the underlying matrix operator.

Any changes to the values of a matrix must be accompanied by a call to NumericFactorization() before the next call to Solve() or the behavior of Solve() is undefined. Any changes to the numerical structure of the matrix must be followed by a call to SymbolicFactorization() and NumericalFactorization() before the next call to Solve().

Once SymbolicFactorization() has been called, classes implementing this interface may assume that any change made to the non-zero structure of the underlying matrix will be accompanied by a call to SymbolicFactorization() prior to a subsequent call to NumericFactorization or Solve().

Named Parameters

Parameters can be changed or added at any time by calling SetParameters(ParamList) with the new parameters specified in ParamList.

It is left to the user to be sure that changes made to the parameters are appropriate for the concrete class that they are using.

Examples of appropriate changes in parameters include:

Examples of inappropriate changes in parameters include:

Results of making inappropriate changes in parameters is unpredictable and could include an error return, a bogus result or ignoring the parameter change.

Transpose solve

Any class implementing Amesos_BaseSolver should handle calls to SetUseTranspose() at any point. However, the result of a call to SetUseTranspose() which is not followed by a call to SymbolicFactorization() and NumericFactorization() is implementation dependent. Some third party libraries are able to solve AT x = b and Ax = b using the same factorization. Others will require a new factorization anytime that a call to SetUseTranspose() changes the intended solve from AT x = b to Ax = b or vice-versa.

Performance expectations

The following is a list of performance guidelines that classes which implement the Amesos_BaseSolver class are expected to maintain.

Memory usage:

For serial codes, no more than one extra copy of the original matrix should be required. Except that some codes require matrix transpostion which requires additional copies of the input matrix.

For distributed memory codes, no serial copies of the original matrix should be required.

Robustness requirements

Failures should be caught by AMESOS_CHK_ERR(). The following error codes should be used:

Because we do not check to see if a matrix has changed between the call to SymbolicFactorization() and the call to NumericFactorization(), it is possible that a change to the matrix will cause a potentially catastrophic error.

Date
Last updated on 24-May-05.

Definition at line 223 of file Amesos_BaseSolver.h.

Constructor & Destructor Documentation

virtual Amesos_BaseSolver::~Amesos_BaseSolver ( )
inlinevirtual

Destructor.

Definition at line 239 of file Amesos_BaseSolver.h.

Member Function Documentation

virtual int Amesos_BaseSolver::SymbolicFactorization ( )
pure 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.

Implemented in Amesos_Merikos, Amesos_Scalapack, Amesos_Klu, Amesos_Mumps, Amesos_Paraklete, Amesos_Superlu, Amesos_Btf, Amesos_Umfpack, Amesos_Lapack, Amesos_Taucs, Amesos_Superludist, Amesos_Dscpack, and Amesos_Pardiso.

virtual int Amesos_BaseSolver::NumericFactorization ( )
pure 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.

Implemented in Amesos_Merikos, Amesos_Scalapack, Amesos_Klu, Amesos_Mumps, Amesos_Paraklete, Amesos_Superlu, Amesos_Btf, Amesos_Umfpack, Amesos_Lapack, Amesos_Taucs, Amesos_Superludist, Amesos_Dscpack, and Amesos_Pardiso.

virtual int Amesos_BaseSolver::Solve ( )
pure 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.

Implemented in Amesos_Merikos, Amesos_Scalapack, Amesos_Klu, Amesos_Mumps, Amesos_Paraklete, Amesos_Btf, Amesos_Superlu, Amesos_Umfpack, Amesos_Lapack, Amesos_Taucs, Amesos_Superludist, Amesos_Dscpack, and Amesos_Pardiso.

virtual int Amesos_BaseSolver::SetUseTranspose ( bool  UseTranspose)
pure virtual

If set true, X will be set to the solution of AT X = B (not A X = B)

If the implementation of this interface does not support transpose use, this method should return a value of -1.

<br >Preconditions:

<br >Postconditions:

  • The next factorization and solve will be performed with the new value of UseTranspose.
Parameters
UseTranspose– (In) If true, solve AT X = B, otherwise solve A X = B.
Returns
Integer error code, set to 0 if successful. Set to -1 if this implementation does not support transpose.

Implemented in Amesos_Merikos, Amesos_Scalapack, Amesos_Klu, Amesos_Btf, Amesos_Mumps, Amesos_Paraklete, Amesos_Superlu, Amesos_Umfpack, Amesos_Taucs, Amesos_Lapack, Amesos_Dscpack, Amesos_Superludist, and Amesos_Pardiso.

virtual bool Amesos_BaseSolver::UseTranspose ( ) const
pure virtual
virtual int Amesos_BaseSolver::SetParameters ( Teuchos::ParameterList ParameterList)
pure 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.

Implemented in Amesos_Merikos, Amesos_Scalapack, Amesos_Btf, Amesos_Klu, Amesos_Mumps, Amesos_Paraklete, Amesos_Superlu, Amesos_Umfpack, Amesos_Lapack, Amesos_Superludist, Amesos_Taucs, Amesos_Dscpack, and Amesos_Pardiso.

virtual const Epetra_LinearProblem* Amesos_BaseSolver::GetProblem ( ) const
pure virtual

Returns the Epetra_LinearProblem.

Warning! Do not call return->SetOperator(...) to attempt to change the Epetra_Operator object (even if the new matrix has the same structure). This new operator matrix will be ignored!

Implemented in Amesos_Mumps, Amesos_Merikos, Amesos_Scalapack, Amesos_Klu, Amesos_Btf, Amesos_Paraklete, Amesos_Superlu, Amesos_Superludist, Amesos_Umfpack, Amesos_Lapack, Amesos_Taucs, Amesos_Dscpack, and Amesos_Pardiso.

virtual bool Amesos_BaseSolver::MatrixShapeOK ( ) const
pure virtual

Returns true if the solver can handle this matrix shape.

Returns true if the matrix shape is one that the underlying sparse direct solver can handle. Classes that work only on square matrices should return false for rectangular matrices. Classes that work only on symmetric matrices whould return false for non-symmetric matrices.

Implemented in Amesos_Mumps, Amesos_Merikos, Amesos_Scalapack, Amesos_Klu, Amesos_Btf, Amesos_Paraklete, Amesos_Superlu, Amesos_Superludist, Amesos_Umfpack, Amesos_Taucs, Amesos_Lapack, Amesos_Dscpack, and Amesos_Pardiso.

virtual const Epetra_Comm& Amesos_BaseSolver::Comm ( ) const
pure virtual
virtual int Amesos_BaseSolver::NumSymbolicFact ( ) const
pure virtual

Returns the number of symbolic factorizations performed by this object.

Implemented in Amesos_Merikos, Amesos_Scalapack, Amesos_Btf, Amesos_Klu, Amesos_Mumps, Amesos_Paraklete, Amesos_Lapack, Amesos_Superlu, Amesos_Umfpack, Amesos_Superludist, Amesos_Taucs, Amesos_Dscpack, and Amesos_Pardiso.

virtual int Amesos_BaseSolver::NumNumericFact ( ) const
pure virtual
virtual int Amesos_BaseSolver::NumSolve ( ) const
pure virtual
virtual void Amesos_BaseSolver::PrintStatus ( ) const
pure virtual
virtual void Amesos_BaseSolver::PrintTiming ( ) const
pure virtual
virtual void Amesos_BaseSolver::setParameterList ( Teuchos::RCP< Teuchos::ParameterList > const &  paramList)
inlinevirtual

Redefined from Teuchos::ParameterListAcceptor (Does Not Work)

Implements Teuchos::ParameterListAcceptor.

Reimplemented in Amesos_Lapack.

Definition at line 409 of file Amesos_BaseSolver.h.

virtual Teuchos::RCP<Teuchos::ParameterList> Amesos_BaseSolver::getNonconstParameterList ( )
inlinevirtual

This is an empty stub.

Implements Teuchos::ParameterListAcceptor.

Definition at line 417 of file Amesos_BaseSolver.h.

virtual Teuchos::RCP<Teuchos::ParameterList> Amesos_BaseSolver::unsetParameterList ( )
inlinevirtual

This is an empty stub.

Implements Teuchos::ParameterListAcceptor.

Reimplemented in Amesos_Lapack.

Definition at line 424 of file Amesos_BaseSolver.h.

virtual void Amesos_BaseSolver::GetTiming ( Teuchos::ParameterList TimingParameterList) const
inlinevirtual

Extracts timing information from the current solver and places it in the parameter list. (Does Not Work)

Reimplemented in Amesos_Scalapack, Amesos_Mumps, Amesos_Klu, Amesos_Paraklete, Amesos_Lapack, Amesos_Superlu, Amesos_Umfpack, Amesos_Superludist, Amesos_Taucs, Amesos_Dscpack, and Amesos_Pardiso.

Definition at line 432 of file Amesos_BaseSolver.h.


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