Amesos
Development
|
Class Amesos_Umfpack: An object-oriented wrapper for UMFPACK. More...
#include <Amesos_Umfpack.h>
Public Member Functions | |
Amesos_Umfpack (const Epetra_LinearProblem &LinearProblem) | |
Amesos_Umfpack Constructor. More... | |
~Amesos_Umfpack (void) | |
Amesos_Umfpack Destructor. More... | |
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_LinearProblem * | GetProblem () const |
Returns the Epetra_LinearProblem. More... | |
bool | MatrixShapeOK () const |
Returns true if UMFPACK can handle this matrix shape. More... | |
int | SetUseTranspose (bool UseTranspose_in) |
If set true, X will be set to the solution of AT X = B (not A X = B) More... | |
bool | UseTranspose () const |
Returns the current UseTranspose setting. | |
const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this operator. | |
double | GetRcond () const |
Returns an estimate of the reciprocal of the condition number. | |
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 from the current solver and places it in the parameter list. | |
Public Member Functions inherited from Amesos_BaseSolver | |
virtual | ~Amesos_BaseSolver () |
Destructor. | |
virtual void | setParameterList (Teuchos::RCP< Teuchos::ParameterList > const ¶mList) |
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 ParameterList > | getParameterList () const |
virtual RCP< const ParameterList > | getValidParameters () const |
Class Amesos_Umfpack: An object-oriented wrapper for UMFPACK.
Amesos_Umfpack will solve a linear systems of equations: A X = B
using Epetra objects and the UMFPACK solver library, where A
is an Epetra_RowMatrix and X
and B
are Epetra_MultiVector objects.
Amesos_Umfpack::Amesos_Umfpack | ( | const Epetra_LinearProblem & | LinearProblem | ) |
Amesos_Umfpack Constructor.
Creates an Amesos_Umfpack instance, using an Epetra_LinearProblem, passing in an already-defined Epetra_LinearProblem object.
Note: The operator in LinearProblem must be an Epetra_RowMatrix.
References SetParameters().
Amesos_Umfpack::~Amesos_Umfpack | ( | void | ) |
Amesos_Umfpack Destructor.
Completely deletes an Amesos_Umfpack object.
References PrintStatus(), Amesos_Status::PrintStatus_, PrintTiming(), Amesos_Status::PrintTiming_, and Amesos_Status::verbose_.
|
inlinevirtual |
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!
Implements Amesos_BaseSolver.
Referenced by Comm(), and MatrixShapeOK().
|
virtual |
Returns true if UMFPACK can handle this matrix shape.
Returns true if the matrix shape is one that UMFPACK can handle. UMFPACK only works with square matrices.
Implements Amesos_BaseSolver.
References GetProblem().
|
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:
<br >Postconditions:
Implements Amesos_BaseSolver.
References Comm(), Amesos_Time::CreateTimer(), Amesos_Status::IsNumericFactorizationOK_, Amesos_Status::IsSymbolicFactorizationOK_, Epetra_Comm::MyPID(), Amesos_Status::NumNumericFact_, Epetra_Comm::NumProc(), and Amesos_Status::NumSymbolicFact_.
Referenced by Solve().
|
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.
Implements Amesos_BaseSolver.
Referenced by Amesos_Umfpack().
|
inlinevirtual |
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:
UseTranspose | – (In) If true, solve AT X = B, otherwise solve A X = B. |
Implements Amesos_BaseSolver.
|
virtual |
Solves A X = B (or AT x = B)
<br >Preconditions:
<br >Postconditions:
Implements Amesos_BaseSolver.
References Amesos_Time::AddTime(), Amesos_Utils::ComputeTrueResidual(), Amesos_Status::ComputeTrueResidual_, Amesos_Utils::ComputeVectorNorms(), Amesos_Status::ComputeVectorNorms_, Epetra_LinearProblem::GetLHS(), Epetra_LinearProblem::GetOperator(), Epetra_LinearProblem::GetRHS(), Insert, Amesos_Status::IsNumericFactorizationOK_, NumericFactorization(), Amesos_Status::NumSolve_, Amesos_Time::ResetTimer(), and UseTranspose().
|
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:
Implements Amesos_BaseSolver.
References Comm(), Amesos_Time::CreateTimer(), Amesos_Status::IsNumericFactorizationOK_, Amesos_Status::IsSymbolicFactorizationOK_, Epetra_Comm::MyPID(), Epetra_Comm::NumProc(), and Amesos_Status::NumSymbolicFact_.