| Amesos Package Browser (Single Doxygen Collection)
    Development
    | 
Amesos_Klu: A serial, unblocked code ideal for getting started and for very sparse matrices, such as circuit matrces. More...
#include <Amesos_Klu.h>

| Private Attributes | |
| int | SerialXlda_ | 
| Teuchos::RCP< Amesos_Klu_Pimpl > | PrivateKluData_ | 
| Teuchos::RCP < Amesos_StandardIndex > | StdIndex_ | 
| Teuchos::RCP < Amesos_StandardIndex > | StdIndexRange_ | 
| Teuchos::RCP < Amesos_StandardIndex > | StdIndexDomain_ | 
| std::vector< int > | Ap | 
| Ap, Ai, Aval form the compressed row storage used by Klu Ai and Aval can point directly into a matrix if it is StorageOptimized(), hence they may either be in vector form or may be a pointer into Epetra_CrsMatrix internals.  More... | |
| std::vector< int > | VecAi | 
| std::vector< double > | VecAval | 
| double * | Aval | 
| int * | Ai | 
| int | UseDataInPlace_ | 
| 1 if Problem_->GetOperator() is stored entirely on process 0  More... | |
| long long | numentries_ | 
| Number of non-zero entries in Problem_->GetOperator()  More... | |
| long long | NumGlobalElements_ | 
| Number of rows and columns in the Problem_->GetOperator()  More... | |
| Epetra_RowMatrix * | RowMatrixA_ | 
| Operator converted to a RowMatrix.  More... | |
| Epetra_CrsMatrix * | CrsMatrixA_ | 
| Operator converted to a CrsMatrix.  More... | |
| Teuchos::RCP< Epetra_Map > | SerialMap_ | 
| Points to a Serial Map (unused if UseDataInPlace_ == 1 )  More... | |
| Teuchos::RCP< Epetra_CrsMatrix > | SerialCrsMatrixA_ | 
| Points to a Serial Copy of A (unused if UseDataInPlace_==1)  More... | |
| Epetra_RowMatrix * | StdIndexMatrix_ | 
| Points to a Contiguous Copy of A.  More... | |
| Epetra_RowMatrix * | SerialMatrix_ | 
| Points to a Serial Copy of A.  More... | |
| bool | TrustMe_ | 
| If true, no checks are made and the matrix is assume to be distributed.  More... | |
| int | NumVectors_ | 
| Number of vectors in RHS and LHS.  More... | |
| double * | SerialXBvalues_ | 
| Pointer to the actual values in the serial version of X and B.  More... | |
| double * | SerialBvalues_ | 
| Teuchos::RCP< Epetra_MultiVector > | SerialB_ | 
| Serial versions of the LHS and RHS (may point to the original vector if serial)  More... | |
| Teuchos::RCP< Epetra_MultiVector > | SerialX_ | 
| Teuchos::RCP< Epetra_MultiVector > | SerialXextract_ | 
| Serial versions of the LHS and RHS (if necessary)  More... | |
| Teuchos::RCP< Epetra_MultiVector > | SerialBextract_ | 
| bool | UseTranspose_ | 
| If true, the transpose of A is used.  More... | |
| const Epetra_LinearProblem * | Problem_ | 
| Pointer to the linear system problem.  More... | |
| std::vector< int > | ColIndicesV_ | 
| Only used for RowMatrices to extract copies.  More... | |
| std::vector< double > | RowValuesV_ | 
| Only used for RowMatrices to extract copies.  More... | |
| Teuchos::RCP< Epetra_Import > | ImportToSerial_ | 
| Importer to process 0.  More... | |
| Teuchos::RCP< Epetra_Import > | ImportRangeToSerial_ | 
| Teuchos::RCP< Epetra_Import > | ImportDomainToSerial_ | 
| int | MtxRedistTime_ | 
| Quick access ids for the individual timings.  More... | |
| int | MtxConvTime_ | 
| int | VecRedistTime_ | 
| int | SymFactTime_ | 
| int | NumFactTime_ | 
| int | SolveTime_ | 
| int | OverheadTime_ | 
|  Private Attributes inherited from Amesos_Control | |
| double | AddToDiag_ | 
| Add thisvalue to the diagonal.  More... | |
| bool | refactorize_ | 
| double | rcond_threshold_ | 
| If error is greater than thisvalue, perform symbolic and numeric factorization with full partial pivoting.  More... | |
| int | ScaleMethod_ | 
| bool | AddZeroToDiag_ | 
| Adds zero to diagonal of redistributed matrix (some solvers choke on a matrix with a partly empty diag)  More... | |
| int | MatrixProperty_ | 
| Set the matrix property.  More... | |
| int | MaxProcesses_ | 
| bool | Reindex_ | 
| If true, the Amesos class should reindex the matrix to standard indexing (i.e.  More... | |
|  Private Attributes inherited from Amesos_Status | |
| bool | IsSymbolicFactorizationOK_ | 
| If true, SymbolicFactorization() has been successfully called.  More... | |
| bool | IsNumericFactorizationOK_ | 
| If true, NumericFactorization() has been successfully called.  More... | |
| bool | PrintTiming_ | 
| If true, prints timing information in the destructor.  More... | |
| bool | PrintStatus_ | 
| If true, print additional information in the destructor.  More... | |
| bool | ComputeVectorNorms_ | 
| If true, prints the norms of X and B in Solve().  More... | |
| bool | ComputeTrueResidual_ | 
| If true, computes the true residual in Solve().  More... | |
| int | verbose_ | 
| Toggles the output level.  More... | |
| int | debug_ | 
| Sets the level of debug_ output.  More... | |
| int | NumSymbolicFact_ | 
| Number of symbolic factorization phases.  More... | |
| int | NumNumericFact_ | 
| Number of numeric factorization phases.  More... | |
| int | NumSolve_ | 
| Number of solves.  More... | |
| double | Threshold_ | 
| int | MyPID_ | 
| int | NumProcs_ | 
| Amesos_Klu (const Epetra_LinearProblem &LinearProblem) | |
| Amesos_Klu Constructor.  More... | |
| ~Amesos_Klu (void) | |
| Amesos_Klu 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 | 
| Get a pointer to the Problem.  More... | |
| 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.  More... | |
| const Epetra_Comm & | Comm () const | 
| Returns a pointer to the Epetra_Comm communicator associated with this operator.  More... | |
| int | SetParameters (Teuchos::ParameterList &ParameterList) | 
| Updates internal variables.  More... | |
| int | NumSymbolicFact () const | 
| Returns the number of symbolic factorizations performed by this object.  More... | |
| int | NumNumericFact () const | 
| Returns the number of numeric factorizations performed by this object.  More... | |
| int | NumSolve () const | 
| Returns the number of solves performed by this object.  More... | |
| void | PrintTiming () const | 
| Prints timing information.  More... | |
| void | PrintStatus () const | 
| Prints information about the factorization and solution phases.  More... | |
| void | GetTiming (Teuchos::ParameterList &TimingParameterList) const | 
| Extracts timing information and places in parameter list.  More... | |
| int | CreateLocalMatrixAndExporters () | 
| int | ExportToSerial () | 
| int | ConvertToKluCRS (bool firsttime) | 
| int | PerformSymbolicFactorization () | 
| int | PerformNumericFactorization () | 
| Additional Inherited Members | |
|  Public Member Functions inherited from Amesos_BaseSolver | |
| virtual | ~Amesos_BaseSolver () | 
| Destructor.  More... | |
| virtual void | setParameterList (Teuchos::RCP< Teuchos::ParameterList > const ¶mList) | 
| 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... | |
|  Public Member Functions inherited from Teuchos::ParameterListAcceptor | |
| virtual RCP< const ParameterList > | getParameterList () const | 
| virtual RCP< const ParameterList > | getValidParameters () const | 
|  Private Member Functions inherited from Amesos_Time | |
| Amesos_Time () | |
| Default constructor to create sizetimers.  More... | |
| virtual | ~Amesos_Time () | 
| Default destructor.  More... | |
| void | CreateTimer (const Epetra_Comm &Comm, int size=1) | 
| Initializes the Time object.  More... | |
| void | ResetTimer (const int timerID=0) | 
| Resets the internally stored time object.  More... | |
| int | AddTime (const std::string what, int dataID, const int timerID=0) | 
| Adds to field whatthe time elapsed since last call to ResetTimer().  More... | |
| double | GetTime (const std::string what) const | 
| Gets the cumulative time using the string.  More... | |
| double | GetTime (const int dataID) const | 
| Gets the cumulative time using the dataID.  More... | |
| void | GetTiming (Teuchos::ParameterList &list) const | 
| Load up the current timing information into the parameter list.  More... | |
|  Private Member Functions inherited from Amesos_NoCopiable | |
| Amesos_NoCopiable () | |
| Default constructor.  More... | |
| ~Amesos_NoCopiable () | |
| Default destructor.  More... | |
|  Private Member Functions inherited from Amesos_Utils | |
| Amesos_Utils () | |
| Default constructor.  More... | |
| ~Amesos_Utils () | |
| Default destructor.  More... | |
| void | ComputeTrueResidual (const Epetra_RowMatrix &Matrix, const Epetra_MultiVector &X, const Epetra_MultiVector &B, const bool UseTranspose, const std::string prefix) const | 
| Computes the true residual, B - Matrix * X, and prints the results.  More... | |
| void | ComputeVectorNorms (const Epetra_MultiVector &X, const Epetra_MultiVector &B, const std::string prefix) const | 
| Computes the norms of X and B and print the results.  More... | |
| void | PrintLine () const | 
| Prints line on std::cout.  More... | |
| void | SetMaxProcesses (int &MaxProcesses, const Epetra_RowMatrix &A) | 
|  Private Member Functions inherited from Amesos_Control | |
| Amesos_Control () | |
| Default constructor.  More... | |
| ~Amesos_Control () | |
| Default destructor.  More... | |
| void | SetControlParameters (const Teuchos::ParameterList &ParameterList) | 
|  Private Member Functions inherited from Amesos_Status | |
| Amesos_Status () | |
| Default constructor.  More... | |
| ~Amesos_Status () | |
| Default destructor.  More... | |
| void | SetStatusParameters (const Teuchos::ParameterList &ParameterList) | 
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:  , where
, where A is an Epetra_RowMatrix and X and B are Epetra_MultiVector objects.
Amesos_Klu computes  more efficiently than
 more efficiently than  . The latter requires a matrix transpose – which costs both time and space.
. 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.
Definition at line 117 of file Amesos_Klu.h.
| 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.
Definition at line 89 of file Amesos_Klu.cpp.
| Amesos_Klu::~Amesos_Klu | ( | void | ) | 
Amesos_Klu Destructor.
Definition at line 111 of file Amesos_Klu.cpp.
| 
 | 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.
Definition at line 617 of file Amesos_Klu.cpp.
| 
 | 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.
Definition at line 681 of file Amesos_Klu.cpp.
| 
 | virtual | 
<br >Preconditions:
<br >Postconditions:
Implements Amesos_BaseSolver.
Definition at line 729 of file Amesos_Klu.cpp.
| 
 | inlinevirtual | 
Get a pointer to the Problem.
Implements Amesos_BaseSolver.
Definition at line 153 of file Amesos_Klu.h.
| 
 | 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.
Definition at line 603 of file Amesos_Klu.cpp.
| 
 | inlinevirtual | 
SetUseTranpose(true) is more efficient in Amesos_Klu.
If SetUseTranspose() is set to true,  is computed.
 is computed. 
Implements Amesos_BaseSolver.
Definition at line 166 of file Amesos_Klu.h.
| 
 | inlinevirtual | 
Returns the current UseTranspose setting.
Implements Amesos_BaseSolver.
Definition at line 168 of file Amesos_Klu.h.
| 
 | inlinevirtual | 
Returns a pointer to the Epetra_Comm communicator associated with this operator.
Implements Amesos_BaseSolver.
Definition at line 170 of file Amesos_Klu.h.
| 
 | 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.
Definition at line 449 of file Amesos_Klu.cpp.
| 
 | inlinevirtual | 
Returns the number of symbolic factorizations performed by this object.
Implements Amesos_BaseSolver.
Definition at line 175 of file Amesos_Klu.h.
| 
 | inlinevirtual | 
Returns the number of numeric factorizations performed by this object.
Implements Amesos_BaseSolver.
Definition at line 178 of file Amesos_Klu.h.
| 
 | inlinevirtual | 
Returns the number of solves performed by this object.
Implements Amesos_BaseSolver.
Definition at line 181 of file Amesos_Klu.h.
| 
 | virtual | 
Prints timing information.
Implements Amesos_BaseSolver.
Definition at line 949 of file Amesos_Klu.cpp.
| 
 | virtual | 
Prints information about the factorization and solution phases.
Implements Amesos_BaseSolver.
Definition at line 926 of file Amesos_Klu.cpp.
| 
 | inlinevirtual | 
Extracts timing information and places in parameter list.
Reimplemented from Amesos_BaseSolver.
Definition at line 190 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 196 of file Amesos_Klu.cpp.
| 
 | private | 
Definition at line 119 of file Amesos_Klu.cpp.
| 
 | private | 
Definition at line 337 of file Amesos_Klu.cpp.
| 
 | private | 
Definition at line 478 of file Amesos_Klu.cpp.
| 
 | private | 
Definition at line 516 of file Amesos_Klu.cpp.
| 
 | private | 
Definition at line 264 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 273 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 274 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 275 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 276 of file Amesos_Klu.h.
| 
 | private | 
Ap, Ai, Aval form the compressed row storage used by Klu Ai and Aval can point directly into a matrix if it is StorageOptimized(), hence they may either be in vector form or may be a pointer into Epetra_CrsMatrix internals.
Ap must always be constructed.
Definition at line 282 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 283 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 284 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 285 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 286 of file Amesos_Klu.h.
| 
 | private | 
1 if Problem_->GetOperator() is stored entirely on process 0
Definition at line 289 of file Amesos_Klu.h.
| 
 | private | 
Number of non-zero entries in Problem_->GetOperator()
Definition at line 291 of file Amesos_Klu.h.
| 
 | private | 
Number of rows and columns in the Problem_->GetOperator()
Definition at line 293 of file Amesos_Klu.h.
| 
 | private | 
Operator converted to a RowMatrix.
Definition at line 296 of file Amesos_Klu.h.
| 
 | private | 
Operator converted to a CrsMatrix.
Definition at line 298 of file Amesos_Klu.h.
| 
 | private | 
Points to a Serial Map (unused if UseDataInPlace_ == 1 )
Definition at line 308 of file Amesos_Klu.h.
| 
 | private | 
Points to a Serial Copy of A (unused if UseDataInPlace_==1)
Definition at line 310 of file Amesos_Klu.h.
| 
 | private | 
Points to a Contiguous Copy of A.
Definition at line 312 of file Amesos_Klu.h.
| 
 | private | 
Points to a Serial Copy of A.
Definition at line 314 of file Amesos_Klu.h.
| 
 | private | 
If true, no checks are made and the matrix is assume to be distributed. 
Definition at line 320 of file Amesos_Klu.h.
| 
 | private | 
Number of vectors in RHS and LHS.
Definition at line 322 of file Amesos_Klu.h.
| 
 | private | 
Pointer to the actual values in the serial version of X and B.
Definition at line 324 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 325 of file Amesos_Klu.h.
| 
 | private | 
Serial versions of the LHS and RHS (may point to the original vector if serial)
Definition at line 327 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 328 of file Amesos_Klu.h.
| 
 | private | 
Serial versions of the LHS and RHS (if necessary)
Definition at line 330 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 331 of file Amesos_Klu.h.
| 
 | private | 
If true, the transpose of A is used. 
Definition at line 334 of file Amesos_Klu.h.
| 
 | private | 
Pointer to the linear system problem.
Definition at line 336 of file Amesos_Klu.h.
| 
 | private | 
Only used for RowMatrices to extract copies.
Definition at line 339 of file Amesos_Klu.h.
| 
 | private | 
Only used for RowMatrices to extract copies.
Definition at line 341 of file Amesos_Klu.h.
| 
 | private | 
Importer to process 0.
Definition at line 343 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 344 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 345 of file Amesos_Klu.h.
| 
 | private | 
Quick access ids for the individual timings.
Definition at line 347 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 347 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 347 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 348 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 348 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 348 of file Amesos_Klu.h.
| 
 | private | 
Definition at line 348 of file Amesos_Klu.h.
 1.8.5
 1.8.5