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

Amesos_Superlu: Amesos interface to Xioye Li's SuperLU 3.0 serial code. More...

#include <Amesos_Superlu.h>

Inheritance diagram for Amesos_Superlu:
Inheritance graph
[legend]

Private Attributes

SLUDatadata_
 Main structure for SuperLU. More...
 
std::vector< double > berr_
 
std::vector< double > ferr_
 
std::vector< int > perm_r_
 
std::vector< int > perm_c_
 
std::vector< int > etree_
 
std::vector< double > R_
 
std::vector< double > C_
 
char equed_
 
double * DummyArray
 stores the matrix in SuperLU format. More...
 
std::vector< int > Ap_
 stores the matrix in SuperLU format. More...
 
std::vector< int > Ai_
 stores the matrix in SuperLU format. More...
 
std::vector< double > Aval_
 
long long NumGlobalRows_
 Global size of the matrix. More...
 
long long NumGlobalNonzeros_
 Global number of nonzeros in the matrix. More...
 
bool UseTranspose_
 If true, solve the linear system with the transpose of the matrix. More...
 
bool FactorizationOK_
 If true, the factorization has been successfully computed. More...
 
bool FactorizationDone_
 
bool ReuseSymbolic_
 
int iam_
 Process number (i.e. Comm().MyPID() More...
 
int MtxConvTime_
 Quick access pointer to internal timing data. More...
 
int MtxRedistTime_
 
int VecRedistTime_
 
int NumFactTime_
 
int SolveTime_
 
int OverheadTime_
 
Teuchos::RCP< Epetra_MapSerialMap_
 Contains a map with all elements assigned to processor 0. More...
 
Teuchos::RCP< Epetra_CrsMatrixSerialCrsMatrixA_
 Contains a matrix with all rows assigned to processor 0. More...
 
Teuchos::RCP< Epetra_ImportImportToSerial_
 Importer from distributed to SerialMap_. More...
 
Epetra_RowMatrixSerialMatrix_
 For parallel runs, stores the matrix defined on SerialMap_. More...
 
const Epetra_LinearProblemProblem_
 Pointer to the user's defined linear problem. More...
 
Epetra_RowMatrixRowMatrixA_
 Pointer to the linear system matrix. More...
 
- 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_
 
 Amesos_Superlu (const Epetra_LinearProblem &LinearProblem)
 Amesos_Superlu Constructor. More...
 
 ~Amesos_Superlu ()
 Amesos_Superlu 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 useTheTranspose)
 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...
 
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 status information. More...
 
void GetTiming (Teuchos::ParameterList &TimingParameterList) const
 Extracts timing information from the current solver and places it in the parameter list. More...
 
const Epetra_MapSerialMap () const
 Returns a reference to the serial map. More...
 
const Epetra_ImportImportToSerial () const
 Returns a reference to the importer. More...
 
int Factor ()
 Factors the matrix, no previous factorization available. More...
 
int ReFactor ()
 Re-factors the matrix. More...
 
int ConvertToSerial ()
 Sets up the matrix on processor 0. More...
 
int PerformNumericFactorization ()
 PerformNumericFactorization - Call Superlu to perform numeric factorization. More...
 

Additional Inherited Members

- Public Member Functions inherited from Amesos_BaseSolver
virtual ~Amesos_BaseSolver ()
 Destructor. More...
 
virtual void setParameterList (Teuchos::RCP< Teuchos::ParameterList > const &paramList)
 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 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)
 

Detailed Description

Amesos_Superlu: Amesos interface to Xioye Li's SuperLU 3.0 serial code.

Class Amesos_Superlu solves the linear systems of equations A X = B, where A is defined as an Epetra_RowMatrix, and X and B are two Epetra_MultiVector's.

Date
Last updated on 28-Apr-05.

Definition at line 92 of file Amesos_Superlu.h.

Constructor & Destructor Documentation

Amesos_Superlu::Amesos_Superlu ( const Epetra_LinearProblem LinearProblem)

Amesos_Superlu Constructor.

Creates an Amesos_Superlu 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 76 of file Amesos_Superlu.cpp.

Amesos_Superlu::~Amesos_Superlu ( void  )

Amesos_Superlu Destructor.

Definition at line 119 of file Amesos_Superlu.cpp.

Member Function Documentation

int Amesos_Superlu::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 409 of file Amesos_Superlu.cpp.

int Amesos_Superlu::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 417 of file Amesos_Superlu.cpp.

int Amesos_Superlu::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 504 of file Amesos_Superlu.cpp.

const Epetra_LinearProblem* Amesos_Superlu::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 127 of file Amesos_Superlu.h.

bool Amesos_Superlu::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 156 of file Amesos_Superlu.cpp.

int Amesos_Superlu::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 131 of file Amesos_Superlu.h.

bool Amesos_Superlu::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Amesos_BaseSolver.

Definition at line 135 of file Amesos_Superlu.h.

const Epetra_Comm& Amesos_Superlu::Comm ( ) const
inlinevirtual

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

Implements Amesos_BaseSolver.

Definition at line 137 of file Amesos_Superlu.h.

int Amesos_Superlu::SetParameters ( Teuchos::ParameterList ParameterList)
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.

Returns
Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

Definition at line 139 of file Amesos_Superlu.cpp.

int Amesos_Superlu::NumSymbolicFact ( ) const
inlinevirtual

Returns the number of symbolic factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 142 of file Amesos_Superlu.h.

int Amesos_Superlu::NumNumericFact ( ) const
inlinevirtual

Returns the number of numeric factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 145 of file Amesos_Superlu.h.

int Amesos_Superlu::NumSolve ( ) const
inlinevirtual

Returns the number of solves performed by this object.

Implements Amesos_BaseSolver.

Definition at line 148 of file Amesos_Superlu.h.

void Amesos_Superlu::PrintTiming ( ) const
virtual

Prints timing information.

Implements Amesos_BaseSolver.

Definition at line 694 of file Amesos_Superlu.cpp.

void Amesos_Superlu::PrintStatus ( ) const
virtual

Prints status information.

Implements Amesos_BaseSolver.

Definition at line 669 of file Amesos_Superlu.cpp.

void Amesos_Superlu::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_Superlu.h.

const Epetra_Map& Amesos_Superlu::SerialMap ( ) const
inlineprivate

Returns a reference to the serial map.

Definition at line 166 of file Amesos_Superlu.h.

const Epetra_Import& Amesos_Superlu::ImportToSerial ( ) const
inlineprivate

Returns a reference to the importer.

Definition at line 173 of file Amesos_Superlu.h.

int Amesos_Superlu::Factor ( void  )
private

Factors the matrix, no previous factorization available.

Definition at line 267 of file Amesos_Superlu.cpp.

int Amesos_Superlu::ReFactor ( )
private

Re-factors the matrix.

Definition at line 332 of file Amesos_Superlu.cpp.

int Amesos_Superlu::ConvertToSerial ( )
private

Sets up the matrix on processor 0.

Definition at line 166 of file Amesos_Superlu.cpp.

int Amesos_Superlu::PerformNumericFactorization ( )
private

PerformNumericFactorization - Call Superlu to perform numeric factorization.

Member Data Documentation

SLUData* Amesos_Superlu::data_
private

Main structure for SuperLU.

Definition at line 193 of file Amesos_Superlu.h.

std::vector<double> Amesos_Superlu::berr_
private

Definition at line 194 of file Amesos_Superlu.h.

std::vector<double> Amesos_Superlu::ferr_
private

Definition at line 195 of file Amesos_Superlu.h.

std::vector<int> Amesos_Superlu::perm_r_
private

Definition at line 196 of file Amesos_Superlu.h.

std::vector<int> Amesos_Superlu::perm_c_
private

Definition at line 197 of file Amesos_Superlu.h.

std::vector<int> Amesos_Superlu::etree_
private

Definition at line 198 of file Amesos_Superlu.h.

std::vector<double> Amesos_Superlu::R_
private

Definition at line 199 of file Amesos_Superlu.h.

std::vector<double> Amesos_Superlu::C_
private

Definition at line 200 of file Amesos_Superlu.h.

char Amesos_Superlu::equed_
private

Definition at line 201 of file Amesos_Superlu.h.

double* Amesos_Superlu::DummyArray
private

stores the matrix in SuperLU format.

Definition at line 203 of file Amesos_Superlu.h.

std::vector<int> Amesos_Superlu::Ap_
private

stores the matrix in SuperLU format.

Definition at line 206 of file Amesos_Superlu.h.

std::vector<int> Amesos_Superlu::Ai_
private

stores the matrix in SuperLU format.

Definition at line 208 of file Amesos_Superlu.h.

std::vector<double> Amesos_Superlu::Aval_
private

Definition at line 210 of file Amesos_Superlu.h.

long long Amesos_Superlu::NumGlobalRows_
private

Global size of the matrix.

Definition at line 212 of file Amesos_Superlu.h.

long long Amesos_Superlu::NumGlobalNonzeros_
private

Global number of nonzeros in the matrix.

Definition at line 214 of file Amesos_Superlu.h.

bool Amesos_Superlu::UseTranspose_
private

If true, solve the linear system with the transpose of the matrix.

Definition at line 216 of file Amesos_Superlu.h.

bool Amesos_Superlu::FactorizationOK_
private

If true, the factorization has been successfully computed.

Definition at line 218 of file Amesos_Superlu.h.

bool Amesos_Superlu::FactorizationDone_
private

Definition at line 219 of file Amesos_Superlu.h.

bool Amesos_Superlu::ReuseSymbolic_
private

Definition at line 220 of file Amesos_Superlu.h.

int Amesos_Superlu::iam_
private

Process number (i.e. Comm().MyPID()

Definition at line 222 of file Amesos_Superlu.h.

int Amesos_Superlu::MtxConvTime_
private

Quick access pointer to internal timing data.

Definition at line 224 of file Amesos_Superlu.h.

int Amesos_Superlu::MtxRedistTime_
private

Definition at line 224 of file Amesos_Superlu.h.

int Amesos_Superlu::VecRedistTime_
private

Definition at line 224 of file Amesos_Superlu.h.

int Amesos_Superlu::NumFactTime_
private

Definition at line 225 of file Amesos_Superlu.h.

int Amesos_Superlu::SolveTime_
private

Definition at line 225 of file Amesos_Superlu.h.

int Amesos_Superlu::OverheadTime_
private

Definition at line 225 of file Amesos_Superlu.h.

Teuchos::RCP<Epetra_Map> Amesos_Superlu::SerialMap_
private

Contains a map with all elements assigned to processor 0.

Definition at line 227 of file Amesos_Superlu.h.

Teuchos::RCP<Epetra_CrsMatrix> Amesos_Superlu::SerialCrsMatrixA_
private

Contains a matrix with all rows assigned to processor 0.

Definition at line 229 of file Amesos_Superlu.h.

Teuchos::RCP<Epetra_Import> Amesos_Superlu::ImportToSerial_
private

Importer from distributed to SerialMap_.

Definition at line 231 of file Amesos_Superlu.h.

Epetra_RowMatrix* Amesos_Superlu::SerialMatrix_
private

For parallel runs, stores the matrix defined on SerialMap_.

Definition at line 233 of file Amesos_Superlu.h.

const Epetra_LinearProblem* Amesos_Superlu::Problem_
private

Pointer to the user's defined linear problem.

Definition at line 235 of file Amesos_Superlu.h.

Epetra_RowMatrix* Amesos_Superlu::RowMatrixA_
private

Pointer to the linear system matrix.

Definition at line 237 of file Amesos_Superlu.h.


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