Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Protected Member Functions | Protected Attributes | List of all members
Amesos_Lapack Class Reference

Amesos_Lapack: an interface to LAPACK. More...

#include <Amesos_Lapack.h>

Inheritance diagram for Amesos_Lapack:
Inheritance graph
[legend]

Protected Member Functions

const Epetra_RowMatrixMatrix () const
 Returns a pointer to the linear system matrix. More...
 
int NumGlobalRows () const
 Returns the number of global rows, or -1 if Matrix() returns 0. More...
 
long long NumGlobalRows64 () const
 
int NumMyRows () const
 Returns the number of local rows, or -1 if Matrix() returns 0. More...
 
const Epetra_MapSerialMap ()
 Returns a reference to serial map (that with all elements on process 0). More...
 
Epetra_RowMatrixSerialMatrix ()
 Returns a reference to serial matrix (that with all rows on process 0). More...
 
Epetra_CrsMatrixSerialCrsMatrix ()
 
const Epetra_ImportMatrixImporter ()
 Returns a reference to the matrix importer (from row map to serial map). More...
 
const Epetra_ExportRhsExporter ()
 Returns a reference to the rhs exporter (from range map to serial map). More...
 
const Epetra_ImportSolutionImporter ()
 Returns a reference to the solution importer (to domain map from serial map). More...
 
int SolveSerial (Epetra_MultiVector &X, const Epetra_MultiVector &B)
 Solves the linear system, when only one process is used. More...
 
int SolveDistributed (Epetra_MultiVector &X, const Epetra_MultiVector &B)
 Solves the linear system, when more than one process is used. More...
 
int DistributedToSerial ()
 Converts a distributed matrix to serial matrix. More...
 
int SerialToDense ()
 Converts a serial matrix to dense format. More...
 
int DenseToFactored ()
 Factors the matrix using LAPACK. More...
 

Protected Attributes

Teuchos::RCP
< Teuchos::ParameterList
pl_
 
Teuchos::RCP< Epetra_RowMatrixSerialMatrix_
 
Teuchos::RCP< Epetra_CrsMatrixSerialCrsMatrix_
 
Teuchos::RCP< Epetra_MapSerialMap_
 
Teuchos::RCP< Epetra_ImportMatrixImporter_
 
Teuchos::RCP< Epetra_ExportRhsExporter_
 
Teuchos::RCP< Epetra_ImportSolutionImporter_
 
Epetra_SerialDenseMatrix DenseMatrix_
 Dense matrix. More...
 
Epetra_SerialDenseMatrix DenseLHS_
 Dense LHS. More...
 
Epetra_SerialDenseMatrix DenseRHS_
 Dense RHS. More...
 
Epetra_SerialDenseSolver DenseSolver_
 Linear problem for dense matrix and vectors. More...
 
bool UseTranspose_
 If true, the linear system with the transpose will be solved. More...
 
const Epetra_LinearProblemProblem_
 Pointer to the linear problem. More...
 
int MtxRedistTime_
 Quick access ids for the individual timings. More...
 
int MtxConvTime_
 
int VecRedistTime_
 
int SymFactTime_
 
int NumFactTime_
 
int SolveTime_
 
long long NumGlobalRows_
 
long long NumGlobalNonzeros_
 
Teuchos::RCP
< Teuchos::ParameterList
ParameterList_
 
 Amesos_Lapack (const Epetra_LinearProblem &LinearProblem)
 Amesos_Lapack Constructor. More...
 
 ~Amesos_Lapack (void)
 Amesos_Lapack 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_LinearProblemGetProblem () const
 Returns the Epetra_LinearProblem. More...
 
bool MatrixShapeOK () const
 Returns true if the solver 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. More...
 
const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this operator. More...
 
void setParameterList (Teuchos::RCP< Teuchos::ParameterList > const &paramList)
 Use this parameter list to read values from. More...
 
Teuchos::RCP
< Teuchos::ParameterList
unsetParameterList ()
 This is an empty stub. More...
 
int SetParameters (Teuchos::ParameterList &ParameterList)
 Deprecated - Sets parameters. More...
 
int GEEV (Epetra_Vector &Er, Epetra_Vector &Ei)
 Computes the eigenvalues of the linear system matrix using DGEEV. 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
 Print timing information. More...
 
void PrintStatus () const
 Print information about the factorization and solution phases. More...
 
void GetTiming (Teuchos::ParameterList &TimingParameterList) const
 Extracts timing information from the current solver and places it in the parameter list. More...
 

Additional Inherited Members

- Public Member Functions inherited from Amesos_BaseSolver
virtual ~Amesos_BaseSolver ()
 Destructor. More...
 
virtual Teuchos::RCP
< Teuchos::ParameterList
getNonconstParameterList ()
 This is an empty stub. More...
 
- Public Member Functions inherited from Teuchos::ParameterListAcceptor
virtual RCP< const ParameterListgetParameterList () const
 
virtual RCP< const ParameterListgetValidParameters () const
 
- Private Member Functions inherited from Amesos_Time
 Amesos_Time ()
 Default constructor to create size timers. 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 what the 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)
 
- Private Attributes inherited from Amesos_Control
double AddToDiag_
 Add this value to the diagonal. More...
 
bool refactorize_
 
double rcond_threshold_
 If error is greater than this value, 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_
 

Detailed Description

Amesos_Lapack: an interface to LAPACK.

Class Amesos_Lapack enables the solution of the distributed linear system, defined by an Epetra_LinearProblem, using LAPACK.

Amesos_Lapack stores the lineaar system matrix as an Epetra_SerialDensMatrix. The linear problem is an Epetra_SerialDenseProblem. Amesos_Lapack factorizes the matrix using DGETRF().

Date
Last updated on 16-Mar-05.
Author
Marzio Sala, 9214.

Definition at line 65 of file Amesos_Lapack.h.

Constructor & Destructor Documentation

Amesos_Lapack::Amesos_Lapack ( const Epetra_LinearProblem LinearProblem)

Amesos_Lapack Constructor.

Creates an Amesos_Lapack 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 39 of file Amesos_Lapack.cpp.

Amesos_Lapack::~Amesos_Lapack ( void  )

Amesos_Lapack Destructor.

Completely deletes an Amesos_Lapack object.

Definition at line 74 of file Amesos_Lapack.cpp.

Member Function Documentation

int Amesos_Lapack::SymbolicFactorization ( )
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.

Implements Amesos_BaseSolver.

Definition at line 127 of file Amesos_Lapack.cpp.

int Amesos_Lapack::NumericFactorization ( )
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.

Implements Amesos_BaseSolver.

Definition at line 183 of file Amesos_Lapack.cpp.

int Amesos_Lapack::Solve ( )
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.

Implements Amesos_BaseSolver.

Definition at line 206 of file Amesos_Lapack.cpp.

const Epetra_LinearProblem* Amesos_Lapack::GetProblem ( ) const
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.

Definition at line 102 of file Amesos_Lapack.h.

bool Amesos_Lapack::MatrixShapeOK ( ) const
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.

Implements Amesos_BaseSolver.

Definition at line 104 of file Amesos_Lapack.cpp.

int Amesos_Lapack::SetUseTranspose ( bool  UseTranspose)
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:

  • 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.

Implements Amesos_BaseSolver.

Definition at line 106 of file Amesos_Lapack.h.

bool Amesos_Lapack::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Amesos_BaseSolver.

Definition at line 111 of file Amesos_Lapack.h.

const Epetra_Comm& Amesos_Lapack::Comm ( ) const
inlinevirtual

Returns a pointer to the Epetra_Comm communicator associated with this operator.

Implements Amesos_BaseSolver.

Definition at line 113 of file Amesos_Lapack.h.

void Amesos_Lapack::setParameterList ( Teuchos::RCP< Teuchos::ParameterList > const &  paramList)
virtual

Use this parameter list to read values from.

Redefined from Teuchos::ParameterListAcceptor

Reimplemented from Amesos_BaseSolver.

Definition at line 60 of file Amesos_Lapack.cpp.

Teuchos::RCP< Teuchos::ParameterList > Amesos_Lapack::unsetParameterList ( )
virtual

This is an empty stub.

Reimplemented from Amesos_BaseSolver.

Definition at line 53 of file Amesos_Lapack.cpp.

int Amesos_Lapack::SetParameters ( Teuchos::ParameterList ParameterList)
virtual

Deprecated - Sets parameters.

Implements Amesos_BaseSolver.

Definition at line 85 of file Amesos_Lapack.cpp.

int Amesos_Lapack::GEEV ( Epetra_Vector Er,
Epetra_Vector Ei 
)

Computes the eigenvalues of the linear system matrix using DGEEV.

Parameters
Er- (Out) On processor zero only, it will contain the real component of the eigenvalues.
Ei- (Out) On processor zero only, it will contain the imaginary component of the eigenvalues.
Note
Er and Ei must have been allocated so that the local length on processor 0 equals the global size of the matrix.

Definition at line 372 of file Amesos_Lapack.cpp.

int Amesos_Lapack::NumSymbolicFact ( ) const
inlinevirtual

Returns the number of symbolic factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 142 of file Amesos_Lapack.h.

int Amesos_Lapack::NumNumericFact ( ) const
inlinevirtual

Returns the number of numeric factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 145 of file Amesos_Lapack.h.

int Amesos_Lapack::NumSolve ( ) const
inlinevirtual

Returns the number of solves performed by this object.

Implements Amesos_BaseSolver.

Definition at line 148 of file Amesos_Lapack.h.

void Amesos_Lapack::PrintTiming ( ) const
virtual

Print timing information.

Implements Amesos_BaseSolver.

Definition at line 488 of file Amesos_Lapack.cpp.

void Amesos_Lapack::PrintStatus ( ) const
virtual

Print information about the factorization and solution phases.

Implements Amesos_BaseSolver.

Definition at line 460 of file Amesos_Lapack.cpp.

void Amesos_Lapack::GetTiming ( Teuchos::ParameterList TimingParameterList) const
inlinevirtual

Extracts timing information from the current solver and places it in the parameter list.

Reimplemented from Amesos_BaseSolver.

Definition at line 157 of file Amesos_Lapack.h.

const Epetra_RowMatrix* Amesos_Lapack::Matrix ( ) const
inlineprotected

Returns a pointer to the linear system matrix.

Definition at line 164 of file Amesos_Lapack.h.

int Amesos_Lapack::NumGlobalRows ( ) const
inlineprotected

Returns the number of global rows, or -1 if Matrix() returns 0.

Definition at line 171 of file Amesos_Lapack.h.

long long Amesos_Lapack::NumGlobalRows64 ( ) const
inlineprotected

Definition at line 177 of file Amesos_Lapack.h.

int Amesos_Lapack::NumMyRows ( ) const
inlineprotected

Returns the number of local rows, or -1 if Matrix() returns 0.

Definition at line 183 of file Amesos_Lapack.h.

const Epetra_Map& Amesos_Lapack::SerialMap ( )
inlineprotected

Returns a reference to serial map (that with all elements on process 0).

Definition at line 189 of file Amesos_Lapack.h.

Epetra_RowMatrix& Amesos_Lapack::SerialMatrix ( )
inlineprotected

Returns a reference to serial matrix (that with all rows on process 0).

Definition at line 195 of file Amesos_Lapack.h.

Epetra_CrsMatrix& Amesos_Lapack::SerialCrsMatrix ( )
inlineprotected

Definition at line 200 of file Amesos_Lapack.h.

const Epetra_Import& Amesos_Lapack::MatrixImporter ( )
inlineprotected

Returns a reference to the matrix importer (from row map to serial map).

Definition at line 206 of file Amesos_Lapack.h.

const Epetra_Export& Amesos_Lapack::RhsExporter ( )
inlineprotected

Returns a reference to the rhs exporter (from range map to serial map).

Definition at line 212 of file Amesos_Lapack.h.

const Epetra_Import& Amesos_Lapack::SolutionImporter ( )
inlineprotected

Returns a reference to the solution importer (to domain map from serial map).

Definition at line 218 of file Amesos_Lapack.h.

int Amesos_Lapack::SolveSerial ( Epetra_MultiVector X,
const Epetra_MultiVector B 
)
protected

Solves the linear system, when only one process is used.

Definition at line 230 of file Amesos_Lapack.cpp.

int Amesos_Lapack::SolveDistributed ( Epetra_MultiVector X,
const Epetra_MultiVector B 
)
protected

Solves the linear system, when more than one process is used.

Definition at line 259 of file Amesos_Lapack.cpp.

int Amesos_Lapack::DistributedToSerial ( )
protected

Converts a distributed matrix to serial matrix.

Definition at line 352 of file Amesos_Lapack.cpp.

int Amesos_Lapack::SerialToDense ( )
protected

Converts a serial matrix to dense format.

Definition at line 304 of file Amesos_Lapack.cpp.

int Amesos_Lapack::DenseToFactored ( )
protected

Factors the matrix using LAPACK.

Definition at line 445 of file Amesos_Lapack.cpp.

Member Data Documentation

Teuchos::RCP<Teuchos::ParameterList> Amesos_Lapack::pl_
protected

Definition at line 223 of file Amesos_Lapack.h.

Teuchos::RCP<Epetra_RowMatrix> Amesos_Lapack::SerialMatrix_
protected

Definition at line 242 of file Amesos_Lapack.h.

Teuchos::RCP<Epetra_CrsMatrix> Amesos_Lapack::SerialCrsMatrix_
protected

Definition at line 243 of file Amesos_Lapack.h.

Teuchos::RCP<Epetra_Map> Amesos_Lapack::SerialMap_
protected

Definition at line 244 of file Amesos_Lapack.h.

Teuchos::RCP<Epetra_Import> Amesos_Lapack::MatrixImporter_
protected

Definition at line 245 of file Amesos_Lapack.h.

Teuchos::RCP<Epetra_Export> Amesos_Lapack::RhsExporter_
protected

Definition at line 246 of file Amesos_Lapack.h.

Teuchos::RCP<Epetra_Import> Amesos_Lapack::SolutionImporter_
protected

Definition at line 247 of file Amesos_Lapack.h.

Epetra_SerialDenseMatrix Amesos_Lapack::DenseMatrix_
protected

Dense matrix.

Definition at line 250 of file Amesos_Lapack.h.

Epetra_SerialDenseMatrix Amesos_Lapack::DenseLHS_
protected

Dense LHS.

Definition at line 252 of file Amesos_Lapack.h.

Epetra_SerialDenseMatrix Amesos_Lapack::DenseRHS_
protected

Dense RHS.

Definition at line 254 of file Amesos_Lapack.h.

Epetra_SerialDenseSolver Amesos_Lapack::DenseSolver_
protected

Linear problem for dense matrix and vectors.

Definition at line 256 of file Amesos_Lapack.h.

bool Amesos_Lapack::UseTranspose_
protected

If true, the linear system with the transpose will be solved.

Definition at line 259 of file Amesos_Lapack.h.

const Epetra_LinearProblem* Amesos_Lapack::Problem_
protected

Pointer to the linear problem.

Definition at line 261 of file Amesos_Lapack.h.

int Amesos_Lapack::MtxRedistTime_
protected

Quick access ids for the individual timings.

Definition at line 264 of file Amesos_Lapack.h.

int Amesos_Lapack::MtxConvTime_
protected

Definition at line 264 of file Amesos_Lapack.h.

int Amesos_Lapack::VecRedistTime_
protected

Definition at line 264 of file Amesos_Lapack.h.

int Amesos_Lapack::SymFactTime_
protected

Definition at line 264 of file Amesos_Lapack.h.

int Amesos_Lapack::NumFactTime_
protected

Definition at line 264 of file Amesos_Lapack.h.

int Amesos_Lapack::SolveTime_
protected

Definition at line 264 of file Amesos_Lapack.h.

long long Amesos_Lapack::NumGlobalRows_
protected

Definition at line 266 of file Amesos_Lapack.h.

long long Amesos_Lapack::NumGlobalNonzeros_
protected

Definition at line 267 of file Amesos_Lapack.h.

Teuchos::RCP<Teuchos::ParameterList> Amesos_Lapack::ParameterList_
protected

Definition at line 269 of file Amesos_Lapack.h.


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