AT X = B more efficiently than
A X = B. More...
#include <Amesos_Klu.h>
Inheritance diagram for Amesos_Klu:
Public Member Functions | |
Constructor methods | |
| Amesos_Klu (const Epetra_LinearProblem &LinearProblem) | |
| Amesos_Klu Constructor. | |
| ~Amesos_Klu (void) | |
| Amesos_Klu Destructor. | |
Mathematical functions. | |
| int | SymbolicFactorization () |
| Performs SymbolicFactorization on the matrix A. | |
| int | NumericFactorization () |
| Performs NumericFactorization on the matrix A. | |
| int | Solve () |
| Solves A X = B (or AT X = B). | |
Additional methods required to support the Epetra_Operator interface. | |
| const Epetra_LinearProblem * | GetProblem () const |
| Get a pointer to the Problem. | |
| bool | MatrixShapeOK () const |
| Returns true if KLU can handle this matrix shape. | |
| int | SetUseTranspose (bool UseTranspose) |
| SetUseTranpose(true) is more efficient in Amesos_Klu. | |
| bool | UseTranspose () const |
| Returns the current UseTranspose setting. | |
| const Epetra_Comm & | Comm () const |
| Returns a pointer to the Epetra_Comm communicator associated with this matrix. | |
| int | SetParameters (Teuchos::ParameterList &ParameterList) |
| Updates internal variables. | |
| void | PrintTiming () |
| Print timing information. | |
| void | PrintStatus () |
| Print information about the factorization and solution phases. | |
Protected Attributes | |
| int * | Lp |
| int * | Li |
| int * | Up |
| int * | Ui |
| int * | P |
| double * | Lx |
| double * | Ux |
| Amesos_Klu_Pimpl * | PrivateKluData_ |
| vector< int > | Ap |
| vector< int > | Ai |
| vector< double > | Aval |
| int | iam |
| int | IsLocal_ |
| int | numentries_ |
| int | NumGlobalElements_ |
| Epetra_Map * | SerialMap_ |
| Epetra_CrsMatrix * | SerialCrsMatrixA_ |
| Epetra_CrsMatrix * | SerialMatrix_ |
| Epetra_CrsMatrix * | TransposeMatrix_ |
| Epetra_CrsMatrix * | Matrix_ |
| bool | UseTranspose_ |
| const Epetra_LinearProblem * | Problem_ |
| bool | PrintTiming_ |
| bool | PrintStatus_ |
| bool | ComputeVectorNorms_ |
| bool | ComputeTrueResidual_ |
| int | verbose_ |
| int | debug_ |
| double | ConTime_ |
| double | SymTime_ |
| double | NumTime_ |
| double | SolTime_ |
| double | VecTime_ |
| double | MatTime_ |
| int | NumSymbolicFact_ |
| int | NumNumericFact_ |
| int | NumSolve_ |
| Epetra_Time | Time |
AT X = B more efficiently than
A X = B.
Amesos_Klu, an object-oriented wrapper for Klu, will solve a linear systems of equations: A X = B using Epetra objects and the Klu solver library, where A is an Epetra_RowMatrix and X and B are Epetra_MultiVector objects.
<br /><br />
AmesosKlu computes
AT X = B more efficiently than
A X = B. The latter requires a matrix transpose - which costs both time and space.
<br /><br />
klu is 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.
<br /><br />
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.
<br /><br />
Klu execution can be tuned through a variety of parameters. Amesos_Klu.h allows control of these parameters through the following named parameters, ignoring parameters with names that it does not recognize. Where possible, the parameters are common to all direct solvers (although some may ignore them). However, some parameters, in particular tuning parameters, are unique to each solver.
|
|
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. |
|
|
Amesos_Klu Destructor. Completely deletes an Amesos_Klu object. |
|
|
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. |
|
|
Performs NumericFactorization on the matrix A. In addition to performing numeric factorization (and symbolic factorization if necessary) 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(). preconditions:
postconditions:
Implements Amesos_BaseSolver. |
|
|
Updates internal variables. <br >Preconditions:
<br >Postconditions:
Implements Amesos_BaseSolver. |
|
|
SetUseTranpose(true) is more efficient in Amesos_Klu.
Implements Amesos_BaseSolver. |
|
|
Solves A X = B (or AT X = B). preconditions:
postconditions:
Implements Amesos_BaseSolver. |
|
|
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(). preconditions:
postconditions:
Implements Amesos_BaseSolver. |
1.3.9.1