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

Amesos_Superludist: An object-oriented wrapper for Superludist. More...

#include <Amesos_Superludist.h>

Inheritance diagram for Amesos_Superludist:
Inheritance graph
[legend]

Public Member Functions

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
 Print various timig. More...
 
void PrintStatus () const
 Print various information about the parameters used by Superludist. More...
 
void GetTiming (Teuchos::ParameterList &TimingParameterList) const
 Extracts timing information from the current solver and places it in the parameter list. More...
 
- 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

const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this operator. More...
 
const Epetra_ImportImporter () const
 
const Epetra_MapUniformMap () const
 
const Epetra_RowMatrixUniformMatrix () const
 
Epetra_CrsMatrixCrsUniformMatrix ()
 
int RedistributeA ()
 
int ReFactor ()
 
int Factor ()
 
- 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

Teuchos::RCP
< Amesos_Superlu_Pimpl
PrivateSuperluData_
 
const Epetra_LinearProblemProblem_
 
Epetra_RowMatrixRowMatrixA_
 
RCP< Epetra_MapUniformMap_
 
RCP< Epetra_CrsMatrixCrsUniformMatrix_
 
RCP< Epetra_RowMatrixUniformMatrix_
 
Teuchos::RCP< Epetra_ImportImporter_
 
bool ReuseSymbolic_
 Allows FactOption to be used on subsequent calls to pdgssvx from NumericFactorization. More...
 
bool Redistribute_
 redistribute the input matrix prior to calling Superludist More...
 
int GridCreated_
 true if the SuperLU_DIST's grid has been created (and has to be free'd) More...
 
int FactorizationDone_
 
bool FactorizationOK_
 true if NumericFactorization() has been successfully called. More...
 
int NumGlobalRows_
 Global dimension of the matrix. More...
 
std::vector< int > Ap_
 
std::vector< int > Ai_
 
std::vector< double > Aval_
 
int * Global_Columns_
 Contains the global ID of local columns. More...
 
int nprow_
 
int npcol_
 
bool PrintNonzeros_
 
std::string ColPerm_
 
std::string RowPerm_
 
int * perm_c_
 
int * perm_r_
 
std::string IterRefine_
 
bool ReplaceTinyPivot_
 
bool Equil_
 
int MtxConvTime_
 
int MtxRedistTime_
 
int VecRedistTime_
 
int NumFactTime_
 
int SolveTime_
 
int OverheadTime_
 
- 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_Superludist (const Epetra_LinearProblem &LinearProblem)
 Amesos_Superludist Constructor. More...
 
 ~Amesos_Superludist (void)
 Amesos_Superludist 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...
 
int SetUseTranspose (bool UseTranspose)
 Amesos_Superludist does not support transpose at this time. More...
 
const Epetra_LinearProblemGetProblem () const
 Returns the Epetra_LinearProblem. More...
 
bool MatrixShapeOK () const
 Returns true if SUPERLUDIST can handle this matrix shape. More...
 
bool UseTranspose () const
 Always returns false. Superludist doesn't support transpose solve. More...
 

Detailed Description

Amesos_Superludist: An object-oriented wrapper for Superludist.

Amesos_Superludist will solve a linear systems of equations: A X = B using Epetra objects and the Superludist solver library, where A is an Epetra_RowMatrix and X and B are Epetra_MultiVector objects.

Definition at line 64 of file Amesos_Superludist.h.

Constructor & Destructor Documentation

Amesos_Superludist::Amesos_Superludist ( const Epetra_LinearProblem LinearProblem)

Amesos_Superludist Constructor.

Creates an Amesos_Superludist 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 116 of file Amesos_Superludist.cpp.

Amesos_Superludist::~Amesos_Superludist ( void  )

Amesos_Superludist Destructor.

Completely deletes an Amesos_Superludist object.

Definition at line 170 of file Amesos_Superludist.cpp.

Member Function Documentation

int Amesos_Superludist::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 688 of file Amesos_Superludist.cpp.

int Amesos_Superludist::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 696 of file Amesos_Superludist.cpp.

int Amesos_Superludist::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 731 of file Amesos_Superludist.cpp.

int Amesos_Superludist::SetUseTranspose ( bool  UseTranspose)
inlinevirtual

Amesos_Superludist does not support transpose at this time.

returns 0 if UseTranspose is set to false, else 1 (failure)

Implements Amesos_BaseSolver.

Definition at line 105 of file Amesos_Superludist.h.

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

bool Amesos_Superludist::MatrixShapeOK ( ) const
virtual

Returns true if SUPERLUDIST can handle this matrix shape.

Returns true if the matrix shape is one that SUPERLUDIST can handle. SUPERLUDIST only works with square matrices.

Implements Amesos_BaseSolver.

Definition at line 678 of file Amesos_Superludist.cpp.

bool Amesos_Superludist::UseTranspose ( ) const
inlinevirtual

Always returns false. Superludist doesn't support transpose solve.

Implements Amesos_BaseSolver.

Definition at line 119 of file Amesos_Superludist.h.

int Amesos_Superludist::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 190 of file Amesos_Superludist.cpp.

int Amesos_Superludist::NumSymbolicFact ( ) const
inlinevirtual

Returns the number of symbolic factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 125 of file Amesos_Superludist.h.

int Amesos_Superludist::NumNumericFact ( ) const
inlinevirtual

Returns the number of numeric factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 128 of file Amesos_Superludist.h.

int Amesos_Superludist::NumSolve ( ) const
inlinevirtual

Returns the number of solves performed by this object.

Implements Amesos_BaseSolver.

Definition at line 131 of file Amesos_Superludist.h.

void Amesos_Superludist::PrintTiming ( ) const
virtual

Print various timig.

Implements Amesos_BaseSolver.

Definition at line 854 of file Amesos_Superludist.cpp.

void Amesos_Superludist::PrintStatus ( ) const
virtual

Print various information about the parameters used by Superludist.

Implements Amesos_BaseSolver.

Definition at line 819 of file Amesos_Superludist.cpp.

void Amesos_Superludist::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 140 of file Amesos_Superludist.h.

const Epetra_Comm& Amesos_Superludist::Comm ( ) const
inlineprivatevirtual

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

Implements Amesos_BaseSolver.

Definition at line 143 of file Amesos_Superludist.h.

const Epetra_Import& Amesos_Superludist::Importer ( ) const
inlineprivate

Definition at line 145 of file Amesos_Superludist.h.

const Epetra_Map& Amesos_Superludist::UniformMap ( ) const
inlineprivate

Definition at line 150 of file Amesos_Superludist.h.

const Epetra_RowMatrix& Amesos_Superludist::UniformMatrix ( ) const
inlineprivate

Definition at line 155 of file Amesos_Superludist.h.

Epetra_CrsMatrix& Amesos_Superludist::CrsUniformMatrix ( )
inlineprivate

Definition at line 160 of file Amesos_Superludist.h.

int Amesos_Superludist::RedistributeA ( )
private

Definition at line 258 of file Amesos_Superludist.cpp.

int Amesos_Superludist::ReFactor ( )
private

Definition at line 590 of file Amesos_Superludist.cpp.

int Amesos_Superludist::Factor ( void  )
private

Definition at line 319 of file Amesos_Superludist.cpp.

Member Data Documentation

Teuchos::RCP<Amesos_Superlu_Pimpl> Amesos_Superludist::PrivateSuperluData_
private

Definition at line 169 of file Amesos_Superludist.h.

const Epetra_LinearProblem* Amesos_Superludist::Problem_
private

Definition at line 176 of file Amesos_Superludist.h.

Epetra_RowMatrix* Amesos_Superludist::RowMatrixA_
private

Definition at line 177 of file Amesos_Superludist.h.

RCP<Epetra_Map> Amesos_Superludist::UniformMap_
private

Definition at line 179 of file Amesos_Superludist.h.

RCP<Epetra_CrsMatrix> Amesos_Superludist::CrsUniformMatrix_
private

Definition at line 180 of file Amesos_Superludist.h.

RCP<Epetra_RowMatrix> Amesos_Superludist::UniformMatrix_
private

Definition at line 181 of file Amesos_Superludist.h.

Teuchos::RCP<Epetra_Import> Amesos_Superludist::Importer_
private

Definition at line 182 of file Amesos_Superludist.h.

bool Amesos_Superludist::ReuseSymbolic_
private

Allows FactOption to be used on subsequent calls to pdgssvx from NumericFactorization.

Definition at line 185 of file Amesos_Superludist.h.

bool Amesos_Superludist::Redistribute_
private

redistribute the input matrix prior to calling Superludist

Definition at line 187 of file Amesos_Superludist.h.

int Amesos_Superludist::GridCreated_
private

true if the SuperLU_DIST's grid has been created (and has to be free'd)

Definition at line 190 of file Amesos_Superludist.h.

int Amesos_Superludist::FactorizationDone_
private

Definition at line 191 of file Amesos_Superludist.h.

bool Amesos_Superludist::FactorizationOK_
private

true if NumericFactorization() has been successfully called.

Definition at line 193 of file Amesos_Superludist.h.

int Amesos_Superludist::NumGlobalRows_
private

Global dimension of the matrix.

Definition at line 196 of file Amesos_Superludist.h.

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

Definition at line 199 of file Amesos_Superludist.h.

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

Definition at line 200 of file Amesos_Superludist.h.

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

Definition at line 201 of file Amesos_Superludist.h.

int* Amesos_Superludist::Global_Columns_
private

Contains the global ID of local columns.

Definition at line 203 of file Amesos_Superludist.h.

int Amesos_Superludist::nprow_
private

Definition at line 205 of file Amesos_Superludist.h.

int Amesos_Superludist::npcol_
private

Definition at line 206 of file Amesos_Superludist.h.

bool Amesos_Superludist::PrintNonzeros_
private

Definition at line 208 of file Amesos_Superludist.h.

std::string Amesos_Superludist::ColPerm_
private

Definition at line 209 of file Amesos_Superludist.h.

std::string Amesos_Superludist::RowPerm_
private

Definition at line 210 of file Amesos_Superludist.h.

int* Amesos_Superludist::perm_c_
private

Definition at line 211 of file Amesos_Superludist.h.

int* Amesos_Superludist::perm_r_
private

Definition at line 212 of file Amesos_Superludist.h.

std::string Amesos_Superludist::IterRefine_
private

Definition at line 213 of file Amesos_Superludist.h.

bool Amesos_Superludist::ReplaceTinyPivot_
private

Definition at line 214 of file Amesos_Superludist.h.

bool Amesos_Superludist::Equil_
private

Definition at line 215 of file Amesos_Superludist.h.

int Amesos_Superludist::MtxConvTime_
private

Definition at line 217 of file Amesos_Superludist.h.

int Amesos_Superludist::MtxRedistTime_
private

Definition at line 217 of file Amesos_Superludist.h.

int Amesos_Superludist::VecRedistTime_
private

Definition at line 217 of file Amesos_Superludist.h.

int Amesos_Superludist::NumFactTime_
private

Definition at line 218 of file Amesos_Superludist.h.

int Amesos_Superludist::SolveTime_
private

Definition at line 218 of file Amesos_Superludist.h.

int Amesos_Superludist::OverheadTime_
private

Definition at line 218 of file Amesos_Superludist.h.


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