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

Amesos_ComponentBaseSolver: A pure virtual class for direct solvers to be used within Amesos_Merikos to form a parallel direct solver. More...

#include <Amesos_ComponentBaseSolver.h>

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

Public Member Functions

virtual ~Amesos_ComponentBaseSolver ()
 Destructor.
 
virtual int PartialFactorization ()=0
 Performs partial factorization on the matrix A. More...
 
virtual int Lsolve ()=0
 Solves L X = B (or LT x = B) More...
 
*virtual int LsolveStart ()=0
 Solves the triangular part of L X1 = B (or LT x = B) More...
 
virtual int LsolvePart (int begin, int end)=0
 Computes L[begin..end,:] X1. More...
 
virtual int Usolve ()=0
 Solves U X = B (or UT x = B) More...
 
*virtual int UsolveStart ()=0
 Solves the triangular part of U X1 = B (or LT x = B) More...
 
virtual int UsolvePart (int begin, int end)=0
 Computes U[:,begin..end] X1. More...
 
virtual int SetRowPermutation (int *RowPermutation)=0
 Solves U X = B (or UT x = B) More...
 
virtual int SetColumnPermutation (int *ColumnPermutation)=0
 SetColumnPermutation.
 
virtual int SetSubMatrixSize (int SubMatrixSize)=0
 SetSubMatrixSize.
 
virtual int GetRowPermutation (int **RowPermutation)=0
 GetRowPermutation. More...
 
virtual int GetColumnPermutation (int **ColumnPermutation)=0
 GetColumnPermutation. More...
 
virtual int GetSubMatrixSize (int *SubMatrixSize)=0
 GetSubMatrixSize.
 
virtual int GetSchurComplement (Epetra_CrsMatrix *SchurComplement)=0
 GetSchurComplement.
 
- Public Member Functions inherited from Amesos_BaseSolver
virtual ~Amesos_BaseSolver ()
 Destructor.
 
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.
 
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.
 
virtual int NumSymbolicFact () const =0
 Returns the number of symbolic factorizations performed by this object.
 
virtual int NumNumericFact () const =0
 Returns the number of numeric factorizations performed by this object.
 
virtual int NumSolve () const =0
 Returns the number of solves performed by this object.
 
virtual void PrintStatus () const =0
 Prints status information about the current solver.
 
virtual void PrintTiming () const =0
 Prints timing information about the current solver.
 
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.
 
virtual void GetTiming (Teuchos::ParameterList &TimingParameterList) const
 Extracts timing information from the current solver and places it in the parameter list. (Does Not Work)
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptor
virtual RCP< const ParameterListgetParameterList () const
 
virtual RCP< const ParameterListgetValidParameters () const
 

Detailed Description

Amesos_ComponentBaseSolver: A pure virtual class for direct solvers to be used within Amesos_Merikos to form a parallel direct solver.

<p>The Amesos_ComponentBaseSolver interface specifies what Amesos_Merikos needs.

Any Amesos class that implements Amesos_ComponentBaseSolver can be used by Amesos_Merikos to perform partial solves on subblocks of the matrix.

<H1>Member functions added by Amesos_ComponentBaseSolver.</H1> 

<ul>
<li>PartialFactorization()
<ul>
   <li>PartialFactorization performs factors at most the 

first SubMatrixSize_ rows and columns. PartialFactorization delays the factorization of any columns which generate unstable (i.e. too small) pivots. PartialFactorization computes and returns the schur complement. PartialFactorization does not need a symbolic factorization phase. It uses the permutation given by SetRowPermutation. Lsolve performs a raw partial solve, treating the unfactored rows and columns as the identity without row or column permutation. Usolve performs a raw partial solve, treating the unfactored rows and columns as the identity without row or column permutation. SetRowPermutation - sets the row permutation GetRowPermutation - gets the row permutation SetColumnPermutation - sets the column permutation GetColumnPermutation - gets the column permutation SetSubMatrixSize - Sets the maximum number of rows (and columns) to factor. GetSubMatrixSize - Returns the number of rows (and columns) actually factored. SchurComplement - Returns the Schur complement, i.e. L21(SubMatrixSize+1:MatrixSize,1:SubMatrixSize) * U12(1:SubMatrixSize,SubMatrixSize+1:MatrixSize)

<H1>Usage Examples</H1> 

<H2>Basic calling sequence</H2> 
    Epetra_LinearProblem Problem(A,X,B);
    Amesos_SolverName Solver(Problem);
    Solver.PartialFactorization() ; 
      ... Ancestor factorization
    Solver.Lsolve() ; 
      ... Ancestor solves
    Solver.Usolve() ; 
<H2>Preconditions:

Postconditions: Constructor requirements

Every Amesos_SolverName class should accept an Epetra_LinearProblem

Member Function Documentation

virtual int Amesos_ComponentBaseSolver::GetColumnPermutation ( int **  ColumnPermutation)
pure virtual

GetColumnPermutation.

ColumnPermutation reflects any row permutations performed by PartialFactorization(). Note: It is not yet clear whether this row permutation includes the ColumnPermuation upon input or whether it returns only the row permuations performed by the most recent call to PartialFactorization(). In other words, in the absence of pivoting, ColumnPermutation might be identical to that given by SetColumnPermutation() or it might be the identity permutation.

virtual int Amesos_ComponentBaseSolver::GetRowPermutation ( int **  RowPermutation)
pure virtual

GetRowPermutation.

RowPermutation reflects any row permutations performed by PartialFactorization(). Note: It is not yet clear whether this row permutation includes the RowPermuation upon input or whether it returns only the row permuations performed by the most recent call to PartialFactorization(). In other words, in the absence of pivoting, RowPermutation might be identical to that given by SetRowPermutation() or it might be the identity permutation.

virtual int Amesos_ComponentBaseSolver::Lsolve ( )
pure virtual

Solves L X = B (or LT x = B)

Returns
Integer error code, set to 0 if successful.
virtual int Amesos_ComponentBaseSolver::LsolvePart ( int  begin,
int  end 
)
pure virtual

Computes L[begin..end,:] X1.

Returns
Integer error code, set to 0 if successful, -1 if unimplimented.
* virtual int Amesos_ComponentBaseSolver::LsolveStart ( )
pure virtual

Solves the triangular part of L X1 = B (or LT x = B)

Returns
Integer error code, set to 0 if successful, -1 if unimplimented.
virtual int Amesos_ComponentBaseSolver::PartialFactorization ( )
pure virtual

Performs partial factorization on the matrix A.

Partial Factorization perfom

Returns
Integer error code, set to 0 if successful.
virtual int Amesos_ComponentBaseSolver::SetRowPermutation ( int *  RowPermutation)
pure virtual

Solves U X = B (or UT x = B)

Returns
Integer error code, set to 0 if successful.SetRowPermutation
virtual int Amesos_ComponentBaseSolver::Usolve ( )
pure virtual

Solves U X = B (or UT x = B)

Returns
Integer error code, set to 0 if successful.
virtual int Amesos_ComponentBaseSolver::UsolvePart ( int  begin,
int  end 
)
pure virtual

Computes U[:,begin..end] X1.

Returns
Integer error code, set to 0 if successful, -1 if unimplimented.
* virtual int Amesos_ComponentBaseSolver::UsolveStart ( )
pure virtual

Solves the triangular part of U X1 = B (or LT x = B)

Returns
Integer error code, set to 0 if successful, -1 if unimplimented.

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