Amesos
Development
|
Amesos_Scalapack: A serial and parallel dense solver. For now, we implement only the unsymmetric ScaLAPACK solver. More...
#include <Amesos_Scalapack.h>
Public Member Functions | |
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_LinearProblem * | GetProblem () const |
Get a pointer to the Problem. | |
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. | |
const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this matrix. | |
int | SetParameters (Teuchos::ParameterList &ParameterList) |
Updates internal variables. More... | |
int | NumSymbolicFact () const |
Returns the number of symbolic factorizations performed by this object. | |
int | NumNumericFact () const |
Returns the number of numeric factorizations performed by this object. | |
int | NumSolve () const |
Returns the number of solves performed by this object. | |
void | PrintTiming () const |
Print timing information. | |
void | PrintStatus () const |
Print information about the factorization and solution phases. | |
void | GetTiming (Teuchos::ParameterList &TimingParameterList) const |
Extracts timing information from the current solver and places it in the parameter list. | |
Public Member Functions inherited from Amesos_BaseSolver | |
virtual | ~Amesos_BaseSolver () |
Destructor. | |
virtual void | setParameterList (Teuchos::RCP< Teuchos::ParameterList > const ¶mList) |
Redefined from Teuchos::ParameterListAcceptor (Does Not Work) | |
virtual Teuchos::RCP < Teuchos::ParameterList > | getNonconstParameterList () |
This is an empty stub. | |
virtual Teuchos::RCP < Teuchos::ParameterList > | unsetParameterList () |
This is an empty stub. | |
Public Member Functions inherited from Teuchos::ParameterListAcceptor | |
virtual RCP< const ParameterList > | getParameterList () const |
virtual RCP< const ParameterList > | getValidParameters () const |
Protected Attributes | |
int | iam_ |
int | NumGlobalElements_ |
int | NumGlobalNonzeros_ |
int | nprow_ |
int | npcol_ |
int | ictxt_ |
int | m_per_p_ |
int | DescA_ [10] |
Epetra_Map * | ScaLAPACK1DMap_ |
Epetra_CrsMatrix * | ScaLAPACK1DMatrix_ |
Epetra_Map * | VectorMap_ |
std::vector< double > | DenseA_ |
std::vector< int > | Ipiv_ |
int | NumOurRows_ |
int | NumOurColumns_ |
bool | UseTranspose_ |
const Epetra_LinearProblem * | Problem_ |
double | ConTime_ |
double | SymTime_ |
double | NumTime_ |
double | SolTime_ |
double | VecTime_ |
double | MatTime_ |
bool | TwoD_distribution_ |
int | grid_nb_ |
int | mypcol_ |
int | myprow_ |
Epetra_CrsMatrix * | FatOut_ |
int | nb_ |
int | lda_ |
Epetra_Time * | Time_ |
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:
There is sufficient memory to store the entrie dense matrix. 8*n^2/p bytes will be required on each process. -AND- one of the following
The matrix is relatively small and dense. Amesos_Scalapack will solve matrices less than 100 by 100 faster than other Amesos classes unless the matrices are very sparse.
The matrix is relatively dense and many processes are available. If a thousand processes are available, Amesos_Scalapack should be competetive with other sparse direct solvers even for matrices whose L and U factors contain only 5% non-zeros.
The matrix is quite dense. Amesos_Scalapack will be well on any matrix whose L and U factors contain 20% or more non-zeros.
Execution time is less important than robustness. Amesos_Scalapack is among the most robust parallel direct solvers.
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") ;
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.
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.
References SetParameters().
Amesos_Scalapack::~Amesos_Scalapack | ( | void | ) |
Amesos_Scalapack Destructor.
Completely deletes an Amesos_Scalapack object.
References PrintStatus(), Amesos_Status::PrintStatus_, PrintTiming(), Amesos_Status::PrintTiming_, and Amesos_Status::verbose_.
|
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.
References GetProblem().
|
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:
Implements Amesos_BaseSolver.
References Comm(), Amesos_Status::debug_, Epetra_LinearProblem::GetOperator(), Epetra_Comm::MyPID(), Epetra_BlockMap::NumGlobalElements(), Epetra_RowMatrix::NumGlobalNonzeros(), Amesos_Status::NumNumericFact_, and Epetra_RowMatrix::RowMatrixRowMap().
|
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.
Implements Amesos_BaseSolver.
References Amesos_Status::debug_, Teuchos::ParameterList::get(), Teuchos::ParameterList::isParameter(), Teuchos::ParameterList::isSublist(), and Teuchos::ParameterList::sublist().
Referenced by Amesos_Scalapack().
|
inlinevirtual |
SetUseTranpose(true) is more efficient in Amesos_Scalapack.
AT X = B
is computed
A X = B
is computed
Implements Amesos_BaseSolver.
References UseTranspose().
|
virtual |
Solves A X = B (or AT X = B)
preconditions:
postconditions:
Implements Amesos_BaseSolver.
References Epetra_Comm::Broadcast(), Comm(), Amesos_Status::ComputeTrueResidual_, Amesos_Status::ComputeVectorNorms_, Amesos_Status::debug_, Epetra_Time::ElapsedTime(), Epetra_LinearProblem::GetLHS(), Epetra_LinearProblem::GetMatrix(), Epetra_LinearProblem::GetOperator(), Epetra_LinearProblem::GetRHS(), Insert, Epetra_RowMatrix::Multiply(), Epetra_Comm::MyPID(), Amesos_Status::NumSolve_, Epetra_Time::ResetStartTime(), Epetra_RowMatrix::RowMatrixRowMap(), UseTranspose(), and Amesos_Status::verbose_.
|
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.
Implements Amesos_BaseSolver.
References Comm(), Amesos_Status::debug_, and Amesos_Status::NumSymbolicFact_.