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

Amesos_Scalapack: A serial and parallel dense solver. For now, we implement only the unsymmetric ScaLAPACK solver. More...

#include <Amesos_Scalapack.h>

Inheritance diagram for Amesos_Scalapack:
Inheritance graph
[legend]

Protected Attributes

int iam_
 
int NumGlobalElements_
 
int NumGlobalNonzeros_
 
int nprow_
 
int npcol_
 
int ictxt_
 
int m_per_p_
 
int DescA_ [10]
 
Epetra_MapScaLAPACK1DMap_
 
Epetra_CrsMatrixScaLAPACK1DMatrix_
 
Epetra_MapVectorMap_
 
std::vector< double > DenseA_
 
std::vector< int > Ipiv_
 
int NumOurRows_
 
int NumOurColumns_
 
bool UseTranspose_
 
const Epetra_LinearProblemProblem_
 
double ConTime_
 
double SymTime_
 
double NumTime_
 
double SolTime_
 
double VecTime_
 
double MatTime_
 
bool TwoD_distribution_
 
int grid_nb_
 
int mypcol_
 
int myprow_
 
Epetra_CrsMatrixFatOut_
 
int nb_
 
int lda_
 
Epetra_TimeTime_
 

Private Member Functions

int RedistributeA ()
 
int ConvertToScalapack ()
 
int PerformNumericFactorization ()
 
- 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)
 
 Amesos_Scalapack (const Epetra_LinearProblem &LinearProblem)
 Amesos_Scalapack Constructor. More...
 
 ~Amesos_Scalapack (void)
 Amesos_Scalapack 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
 Get a pointer to the Problem. More...
 
bool MatrixShapeOK () const
 Returns true if SCALAPACK can handle this matrix shape. More...
 
int SetUseTranspose (bool UseTranspose)
 SetUseTranpose(true) is more efficient in Amesos_Scalapack. 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 matrix. 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
 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 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 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_Scalapack: A serial and parallel dense solver. For now, we implement only the unsymmetric ScaLAPACK solver.

 Amesos_Scalapack, an object-oriented wrapper for LAPACK and ScaLAPACK, 

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



Amesos_Scalapack can be competitive for matrices that are not particularly sparse. ScaLAPACK solves matrices for which the fill-in is roughly 10% to 20% of the matrix size in time comparable to that achieve by other Amesos classes. Amesos_Scalapack scales well and hence its performance advantage will be largest when large number of processes are involved.



Amesos_Scalapack uses the ScaLAPACK functions PDGETRF and PDGETRS if more than one process is used. If only one process is used, Amesos_ScaLAPACK uses the LAPACK function PDGETRF and PDGETRS.



AmesosScaLAPACK uses full partial pivoting and will therefore provide answers that are at least as accurate as any direct sparse solver.



AmesosScalapack makes sense under the following circumstances:

Common control parameters :

Amesos_Scalapack supports the following parameters which are common to across multiple Amesos solvers:

Amesos_Scalapack supports the following parameters specific to Amesos_Scalapack.
Teuchos::ParameterList ScalapackParams = ParameterList.sublist("Scalapack") ;

Limitations:

None of the following limitations would be particularly difficult to remove.



The present implementation limits the number of right hand sides to the number of rows assigned to each process. i.e. nrhs < n/p.



The present implementation does not take advantage of symmetric or symmetric positive definite matrices, although ScaLAPACK has separate routines to take advantages of such matrices.

Definition at line 146 of file Amesos_Scalapack.h.

Constructor & Destructor Documentation

Amesos_Scalapack::Amesos_Scalapack ( const Epetra_LinearProblem LinearProblem)

Amesos_Scalapack Constructor.

Creates an Amesos_Scalapack 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 50 of file Amesos_Scalapack.cpp.

Amesos_Scalapack::~Amesos_Scalapack ( void  )

Amesos_Scalapack Destructor.

Completely deletes an Amesos_Scalapack object.

Definition at line 71 of file Amesos_Scalapack.cpp.

Member Function Documentation

int Amesos_Scalapack::SymbolicFactorization ( )
virtual

Performs SymbolicFactorization on the matrix A.

There is no symbolic factorization phase in ScaLAPACK, as it operates only on dense matrices. Hence, Amesos_Scalapack::SymbolicFactorization() takes no action.

Returns
Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

Definition at line 778 of file Amesos_Scalapack.cpp.

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

preconditions:

postconditions:

  • nprow_, npcol_, DescA_
  • DenseA will be factored
  • Ipiv_ contains the pivots
Returns
Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

Definition at line 789 of file Amesos_Scalapack.cpp.

int Amesos_Scalapack::Solve ( void  )
virtual

Solves A X = B (or AT X = B)

preconditions:

postconditions:

  • X will be set such that A X = B (or AT X = B), within the limits of the accuracy of the the Scalapack solver.
Returns
Integer error code, set to 0 if successful.

Implements Amesos_BaseSolver.

Definition at line 810 of file Amesos_Scalapack.cpp.

const Epetra_LinearProblem* Amesos_Scalapack::GetProblem ( ) const
inlinevirtual

Get a pointer to the Problem.

Implements Amesos_BaseSolver.

Definition at line 246 of file Amesos_Scalapack.h.

bool Amesos_Scalapack::MatrixShapeOK ( ) const
virtual

Returns true if SCALAPACK can handle this matrix shape.

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

Implements Amesos_BaseSolver.

Definition at line 769 of file Amesos_Scalapack.cpp.

int Amesos_Scalapack::SetUseTranspose ( bool  UseTranspose)
inlinevirtual

SetUseTranpose(true) is more efficient in Amesos_Scalapack.

Implements Amesos_BaseSolver.

Definition at line 267 of file Amesos_Scalapack.h.

bool Amesos_Scalapack::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Amesos_BaseSolver.

Definition at line 270 of file Amesos_Scalapack.h.

const Epetra_Comm& Amesos_Scalapack::Comm ( ) const
inlinevirtual

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

Implements Amesos_BaseSolver.

Definition at line 273 of file Amesos_Scalapack.h.

int Amesos_Scalapack::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 subsequent 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 657 of file Amesos_Scalapack.cpp.

int Amesos_Scalapack::NumSymbolicFact ( ) const
inlinevirtual

Returns the number of symbolic factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 295 of file Amesos_Scalapack.h.

int Amesos_Scalapack::NumNumericFact ( ) const
inlinevirtual

Returns the number of numeric factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 298 of file Amesos_Scalapack.h.

int Amesos_Scalapack::NumSolve ( ) const
inlinevirtual

Returns the number of solves performed by this object.

Implements Amesos_BaseSolver.

Definition at line 301 of file Amesos_Scalapack.h.

void Amesos_Scalapack::PrintTiming ( ) const
virtual

Print timing information.

Implements Amesos_BaseSolver.

Definition at line 1015 of file Amesos_Scalapack.cpp.

void Amesos_Scalapack::PrintStatus ( ) const
virtual

Print information about the factorization and solution phases.

Implements Amesos_BaseSolver.

Definition at line 994 of file Amesos_Scalapack.cpp.

void Amesos_Scalapack::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 310 of file Amesos_Scalapack.h.

int Amesos_Scalapack::RedistributeA ( )
private

Definition at line 124 of file Amesos_Scalapack.cpp.

int Amesos_Scalapack::ConvertToScalapack ( )
private

Definition at line 510 of file Amesos_Scalapack.cpp.

int Amesos_Scalapack::PerformNumericFactorization ( )
private

Definition at line 694 of file Amesos_Scalapack.cpp.

Member Data Documentation

int Amesos_Scalapack::iam_
protected

Definition at line 355 of file Amesos_Scalapack.h.

int Amesos_Scalapack::NumGlobalElements_
protected

Definition at line 357 of file Amesos_Scalapack.h.

int Amesos_Scalapack::NumGlobalNonzeros_
protected

Definition at line 358 of file Amesos_Scalapack.h.

int Amesos_Scalapack::nprow_
protected

Definition at line 363 of file Amesos_Scalapack.h.

int Amesos_Scalapack::npcol_
protected

Definition at line 364 of file Amesos_Scalapack.h.

int Amesos_Scalapack::ictxt_
protected

Definition at line 365 of file Amesos_Scalapack.h.

int Amesos_Scalapack::m_per_p_
protected

Definition at line 366 of file Amesos_Scalapack.h.

int Amesos_Scalapack::DescA_[10]
protected

Definition at line 367 of file Amesos_Scalapack.h.

Epetra_Map* Amesos_Scalapack::ScaLAPACK1DMap_
protected

Definition at line 369 of file Amesos_Scalapack.h.

Epetra_CrsMatrix* Amesos_Scalapack::ScaLAPACK1DMatrix_
protected

Definition at line 371 of file Amesos_Scalapack.h.

Epetra_Map* Amesos_Scalapack::VectorMap_
protected

Definition at line 373 of file Amesos_Scalapack.h.

std::vector<double> Amesos_Scalapack::DenseA_
protected

Definition at line 374 of file Amesos_Scalapack.h.

std::vector<int> Amesos_Scalapack::Ipiv_
protected

Definition at line 375 of file Amesos_Scalapack.h.

int Amesos_Scalapack::NumOurRows_
protected

Definition at line 376 of file Amesos_Scalapack.h.

int Amesos_Scalapack::NumOurColumns_
protected

Definition at line 377 of file Amesos_Scalapack.h.

bool Amesos_Scalapack::UseTranspose_
protected

Definition at line 380 of file Amesos_Scalapack.h.

const Epetra_LinearProblem* Amesos_Scalapack::Problem_
protected

Definition at line 381 of file Amesos_Scalapack.h.

double Amesos_Scalapack::ConTime_
protected

Definition at line 385 of file Amesos_Scalapack.h.

double Amesos_Scalapack::SymTime_
protected

Definition at line 386 of file Amesos_Scalapack.h.

double Amesos_Scalapack::NumTime_
protected

Definition at line 387 of file Amesos_Scalapack.h.

double Amesos_Scalapack::SolTime_
protected

Definition at line 388 of file Amesos_Scalapack.h.

double Amesos_Scalapack::VecTime_
protected

Definition at line 389 of file Amesos_Scalapack.h.

double Amesos_Scalapack::MatTime_
protected

Definition at line 390 of file Amesos_Scalapack.h.

bool Amesos_Scalapack::TwoD_distribution_
protected

Definition at line 395 of file Amesos_Scalapack.h.

int Amesos_Scalapack::grid_nb_
protected

Definition at line 396 of file Amesos_Scalapack.h.

int Amesos_Scalapack::mypcol_
protected

Definition at line 397 of file Amesos_Scalapack.h.

int Amesos_Scalapack::myprow_
protected

Definition at line 398 of file Amesos_Scalapack.h.

Epetra_CrsMatrix* Amesos_Scalapack::FatOut_
protected

Definition at line 399 of file Amesos_Scalapack.h.

int Amesos_Scalapack::nb_
protected

Definition at line 404 of file Amesos_Scalapack.h.

int Amesos_Scalapack::lda_
protected

Definition at line 405 of file Amesos_Scalapack.h.

Epetra_Time* Amesos_Scalapack::Time_
protected

Definition at line 407 of file Amesos_Scalapack.h.


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