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

Amesos_Mumps: An object-oriented wrapper for the double precision version of MUMPS. More...

#include <Amesos_Mumps.h>

Inheritance diagram for Amesos_Mumps:
Inheritance graph
[legend]

Public Member Functions

bool MatrixShapeOK () const
 Returns true if the solver can handle this matrix shape. More...
 
const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this matrix. More...
 
const Epetra_LinearProblemGetProblem () const
 Gets a pointer to the Epetra_LinearProblem. 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
 

Protected Member Functions

Epetra_RowMatrixMatrix ()
 Returns a reference to the linear system matrix. More...
 
const Epetra_RowMatrixMatrix () const
 
Epetra_MapRedistrMap ()
 Returns a reference to the map for redistributed matrix. More...
 
Epetra_ImportRedistrImporter ()
 Returns a reference for the redistributed importer. More...
 
Epetra_RowMatrixRedistrMatrix (const bool ImportMatrix=false)
 Returns a reference to the redistributed matrix, imports it is ImportMatrix is true. More...
 
Epetra_MapSerialMap ()
 Returns a reference to the map with all elements on process 0. More...
 
Epetra_ImportSerialImporter ()
 Returns a reference to the importer for SerialMap(). More...
 
int ConvertToTriplet (const bool OnlyValues)
 Converts to MUMPS format (COO format). More...
 
int CheckError ()
 Checks for MUMPS error, prints them if any. See MUMPS' manual. More...
 
void CheckParameters ()
 Check input parameters. More...
 
void SetICNTLandCNTL ()
 

Protected Attributes

bool IsConvertToTripletOK_
 true if matrix has already been converted to COO format More...
 
bool IsComputeSchurComplementOK_
 true if the Schur complement has been computed (need to free memory) More...
 
bool NoDestroy_
 
std::vector< int > Row
 row indices of nonzero elements More...
 
std::vector< int > Col
 column indices of nonzero elements More...
 
std::vector< double > Val
 values of nonzero elements More...
 
int MaxProcs_
 Maximum number of processors in the MUMPS' communicator. More...
 
bool UseTranspose_
 If true, solve the problem with AT. More...
 
int MtxConvTime_
 Quick access pointers to the internal timers. More...
 
int MtxRedistTime_
 
int VecRedistTime_
 
int SymFactTime_
 
int NumFactTime_
 
int SolveTime_
 
double * RowSca_
 Row and column scaling. More...
 
double * ColSca_
 
int * PermIn_
 PermIn for MUMPS. More...
 
int NumSchurComplementRows_
 Number of rows in the Schur complement (if required) More...
 
int * SchurComplementRows_
 Rows for the Schur complement (if required) More...
 
RCP< Epetra_CrsMatrixCrsSchurComplement_
 Pointer to the Schur complement, as CrsMatrix. More...
 
RCP< Epetra_SerialDenseMatrixDenseSchurComplement_
 Pointer to the Schur complement,as DenseMatrix. More...
 
const Epetra_LinearProblemProblem_
 Pointer to the linear problem to be solved. More...
 
RCP< Epetra_MapRedistrMap_
 Redistributed matrix. More...
 
RCP< Epetra_ImportRedistrImporter_
 Redistributed importer (from Matrix().RowMatrixRowMap() to RedistrMatrix().RowMatrixRowMap()). More...
 
RCP< Epetra_CrsMatrixRedistrMatrix_
 Redistributed matrix (only if MaxProcs_ > 1). More...
 
RCP< Epetra_MapSerialMap_
 Map with all elements on process 0 (for solution and rhs). More...
 
RCP< Epetra_ImportSerialImporter_
 Importer from Matrix.OperatorDomainMap() to SerialMap_. More...
 
DMUMPS_STRUC_C MDS
 
std::map< int, int > ICNTL
 
std::map< int, double > CNTL
 
 Amesos_Mumps (const Epetra_LinearProblem &LinearProblem)
 Amesos_Mumps Constructor. More...
 
 ~Amesos_Mumps (void)
 Amesos_Mumps 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...
 
void Destroy ()
 Destroys all data associated with this object. 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...
 
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 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...
 
int SetPrecscaling (double *ColSca, double *RowSca)
 Set prescaling. More...
 
int SetRowScaling (double *RowSca)
 Set row scaling. More...
 
int SetColScaling (double *ColSca)
 Set column scaling. More...
 
int SetOrdering (int *PermIn)
 Sets ordering. More...
 
double * GetRINFO ()
 Gets the pointer to the RINFO array (defined on all processes). More...
 
int * GetINFO ()
 Gets the pointer to the INFO array (defined on all processes). More...
 
double * GetRINFOG ()
 Gets the pointer to the RINFOG array (defined on host only). More...
 
int * GetINFOG ()
 Get the pointer to the INFOG array (defined on host only). More...
 
void SetICNTL (int pos, int value)
 Set ICNTL[pos] to value. pos is expressed in FORTRAN style (starting from 1). More...
 
void SetCNTL (int pos, double value)
 Set CNTL[pos] to value. pos is expressed in FORTRAN style (starting from 1). More...
 

Additional Inherited Members

- 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_Mumps: An object-oriented wrapper for the double precision version of MUMPS.

Amesos_Mumps is an interface to the the double precision version of the sparse parallel direct solver MUMPS. Given an Epetra_RowMatrix A, and two Epetra_MultiVectors X and B, the solution with Amesos_Mumps reads as follows:

  1. Epetra_LinearProblem Problem; Amesos_BaseSolver * Solver; Amesos Amesos_Factory;
  2. Solver = Amesos_Factory.Create("Amesos_Mumps", Problem);
  3. if( Solver == 0 ) std::cerr << "library not available" << std::endl;
  4. Problem.SetMatrix(&A);
  5. Solver->SymbolicFactorization();
  6. Solver->NumericFactorization();
  7. Problem.SetLHS(&X);
  8. Problem.SetLHS(&B);
  9. Solver->Solve();

A number of parameters is available to tune the performances of MUMPS. We refer to the Amesos Reference Guide for a detailed overview of these parameters. Here, we just recall that it is possible to solve the linear system on a subset of the processes contained in the Comm object of the Epetra_LinearProblem.

Amesos_Mumps accepts any Epetra_RowMatrix derived class. However, special functions are available for Epetra_CrsMatrix and Epetra_VbrMatrix objects.

As Amesos is based on Epetra, and Epetra is only double-precision, we still require an Epetra_LinearProblem composed by a double-precision matrix, and two double-precision vectors. The solution vector is casted to double after solution. Single precision may be of interest if Amesos is used with ML, to solve the coarse problem (for which single-precision can be enough in term of numerical error, and usually save memory and CPU time).

Amesos_Mumps is based on Amesos_EpetraBaseSolver, that is derived from Amesos_BaseSolver. The main redistribution utilities, as well as a getrow function, is obtained by EpetraBaseSolver.

Warning
This interface is compatible with MUMPS 4.5.4.
Date
Last modified 26-Jan-06
Author
Marzio Sala, ETHZ.

Definition at line 112 of file Amesos_Mumps.h.

Constructor & Destructor Documentation

Amesos_Mumps::Amesos_Mumps ( const Epetra_LinearProblem LinearProblem)

Amesos_Mumps Constructor.

Creates an Amesos_Mumps instance, using an Epetra_LinearProblem,

Definition at line 51 of file Amesos_Mumps.cpp.

Amesos_Mumps::~Amesos_Mumps ( void  )

Amesos_Mumps Destructor.

Deletes an Amesos_Mumps object.

Definition at line 165 of file Amesos_Mumps.cpp.

Member Function Documentation

int Amesos_Mumps::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 415 of file Amesos_Mumps.cpp.

int Amesos_Mumps::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 574 of file Amesos_Mumps.cpp.

int Amesos_Mumps::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 620 of file Amesos_Mumps.cpp.

void Amesos_Mumps::Destroy ( )

Destroys all data associated with this object.

Definition at line 134 of file Amesos_Mumps.cpp.

int Amesos_Mumps::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 145 of file Amesos_Mumps.h.

bool Amesos_Mumps::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Amesos_BaseSolver.

Definition at line 153 of file Amesos_Mumps.h.

int Amesos_Mumps::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 299 of file Amesos_Mumps.cpp.

int Amesos_Mumps::NumSymbolicFact ( ) const
inlinevirtual

Returns the number of symbolic factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 162 of file Amesos_Mumps.h.

int Amesos_Mumps::NumNumericFact ( ) const
inlinevirtual

Returns the number of numeric factorizations performed by this object.

Implements Amesos_BaseSolver.

Definition at line 165 of file Amesos_Mumps.h.

int Amesos_Mumps::NumSolve ( ) const
inlinevirtual

Returns the number of solves performed by this object.

Implements Amesos_BaseSolver.

Definition at line 168 of file Amesos_Mumps.h.

void Amesos_Mumps::PrintTiming ( ) const
virtual

Prints timing information.

In the destruction phase, prints out detailed information about the various phases: symbolic and numeric factorization, solution, gather/scatter for vectors and matrices.

Implements Amesos_BaseSolver.

Definition at line 843 of file Amesos_Mumps.cpp.

void Amesos_Mumps::PrintStatus ( ) const
virtual

Prints information about the factorization and solution phases.

In the destruction phase, prints out some information furnished by MUMPS, like the amount of required memory, the MFLOPS.

Implements Amesos_BaseSolver.

Definition at line 767 of file Amesos_Mumps.cpp.

void Amesos_Mumps::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 188 of file Amesos_Mumps.h.

int Amesos_Mumps::SetPrecscaling ( double *  ColSca,
double *  RowSca 
)
inline

Set prescaling.

Use double precision vectors of size N (global dimension of the matrix) as scaling for columns and rows. ColSca and RowSca must be defined on the host only, and allocated by the user, if the user sets ICNTL(8) = -1.

Both input vectors are float with –enable-amesos-smumps, double otherwise.

Definition at line 233 of file Amesos_Mumps.h.

int Amesos_Mumps::SetRowScaling ( double *  RowSca)
inline

Set row scaling.

Use double precision vectors of size N (global dimension of the matrix) for row scaling.

Parameters
RowSca(In) -float pointer with –enable-amesos-smumps, double pointer otherwise.

Definition at line 246 of file Amesos_Mumps.h.

int Amesos_Mumps::SetColScaling ( double *  ColSca)
inline

Set column scaling.

Use double precision vectors of size N (global dimension of the matrix) for column scaling.

Parameters
ColSca(In) - float pointer with –enable-amesos-smumps, double pointer otherwise.

Definition at line 258 of file Amesos_Mumps.h.

int Amesos_Mumps::SetOrdering ( int *  PermIn)
inline

Sets ordering.

Use integer vectors of size N (global dimension of the matrix) as given ordering. PermIn must be defined on the host only, and allocated by the user, if the user sets ICNTL(7) = 1.

Definition at line 269 of file Amesos_Mumps.h.

double * Amesos_Mumps::GetRINFO ( )

Gets the pointer to the RINFO array (defined on all processes).

Gets the pointer to the internally stored RINFO array, of type float if option –enable-amesos-smumps is enabled, double otherwise.

Definition at line 975 of file Amesos_Mumps.cpp.

int * Amesos_Mumps::GetINFO ( )

Gets the pointer to the INFO array (defined on all processes).

Gets the pointer to the internally stored INFO array, of type int.

Definition at line 981 of file Amesos_Mumps.cpp.

double * Amesos_Mumps::GetRINFOG ( )

Gets the pointer to the RINFOG array (defined on host only).

Gets the pointer to the internally stored RINFOG array (defined on the host process only), of type float if option –enable-amesos-smumps is enabled, double otherwise.

Definition at line 987 of file Amesos_Mumps.cpp.

int * Amesos_Mumps::GetINFOG ( )

Get the pointer to the INFOG array (defined on host only).

Gets the pointer to the internally stored INFOG (defined on the host process only) array, of type int.

Definition at line 993 of file Amesos_Mumps.cpp.

void Amesos_Mumps::SetICNTL ( int  pos,
int  value 
)

Set ICNTL[pos] to value. pos is expressed in FORTRAN style (starting from 1).

Definition at line 241 of file Amesos_Mumps.cpp.

void Amesos_Mumps::SetCNTL ( int  pos,
double  value 
)

Set CNTL[pos] to value. pos is expressed in FORTRAN style (starting from 1).

Definition at line 248 of file Amesos_Mumps.cpp.

bool Amesos_Mumps::MatrixShapeOK ( ) const
inlinevirtual

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 308 of file Amesos_Mumps.h.

const Epetra_Comm& Amesos_Mumps::Comm ( ) const
inlinevirtual

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

Implements Amesos_BaseSolver.

Definition at line 319 of file Amesos_Mumps.h.

const Epetra_LinearProblem* Amesos_Mumps::GetProblem ( ) const
inlinevirtual

Gets a pointer to the Epetra_LinearProblem.

Implements Amesos_BaseSolver.

Definition at line 322 of file Amesos_Mumps.h.

Epetra_RowMatrix & Amesos_Mumps::Matrix ( )
protected

Returns a reference to the linear system matrix.

Definition at line 885 of file Amesos_Mumps.cpp.

const Epetra_RowMatrix & Amesos_Mumps::Matrix ( ) const
protected

Definition at line 893 of file Amesos_Mumps.cpp.

Epetra_Map & Amesos_Mumps::RedistrMap ( )
protected

Returns a reference to the map for redistributed matrix.

Definition at line 901 of file Amesos_Mumps.cpp.

Epetra_Import & Amesos_Mumps::RedistrImporter ( )
protected

Returns a reference for the redistributed importer.

Definition at line 919 of file Amesos_Mumps.cpp.

Epetra_RowMatrix & Amesos_Mumps::RedistrMatrix ( const bool  ImportMatrix = false)
protected

Returns a reference to the redistributed matrix, imports it is ImportMatrix is true.

Definition at line 931 of file Amesos_Mumps.cpp.

Epetra_Map & Amesos_Mumps::SerialMap ( )
protected

Returns a reference to the map with all elements on process 0.

Definition at line 952 of file Amesos_Mumps.cpp.

Epetra_Import & Amesos_Mumps::SerialImporter ( )
protected

Returns a reference to the importer for SerialMap().

Definition at line 965 of file Amesos_Mumps.cpp.

int Amesos_Mumps::ConvertToTriplet ( const bool  OnlyValues)
protected

Converts to MUMPS format (COO format).

Definition at line 171 of file Amesos_Mumps.cpp.

int Amesos_Mumps::CheckError ( )
protected

Checks for MUMPS error, prints them if any. See MUMPS' manual.

Definition at line 814 of file Amesos_Mumps.cpp.

void Amesos_Mumps::CheckParameters ( )
protected

Check input parameters.

Definition at line 372 of file Amesos_Mumps.cpp.

void Amesos_Mumps::SetICNTLandCNTL ( )
protected

Definition at line 254 of file Amesos_Mumps.cpp.

Member Data Documentation

bool Amesos_Mumps::IsConvertToTripletOK_
protected

true if matrix has already been converted to COO format

Definition at line 358 of file Amesos_Mumps.h.

bool Amesos_Mumps::IsComputeSchurComplementOK_
protected

true if the Schur complement has been computed (need to free memory)

Definition at line 360 of file Amesos_Mumps.h.

bool Amesos_Mumps::NoDestroy_
protected

Definition at line 363 of file Amesos_Mumps.h.

std::vector<int> Amesos_Mumps::Row
protected

row indices of nonzero elements

Definition at line 366 of file Amesos_Mumps.h.

std::vector<int> Amesos_Mumps::Col
protected

column indices of nonzero elements

Definition at line 368 of file Amesos_Mumps.h.

std::vector<double> Amesos_Mumps::Val
protected

values of nonzero elements

Definition at line 370 of file Amesos_Mumps.h.

int Amesos_Mumps::MaxProcs_
protected

Maximum number of processors in the MUMPS' communicator.

Definition at line 373 of file Amesos_Mumps.h.

bool Amesos_Mumps::UseTranspose_
protected

If true, solve the problem with AT.

Definition at line 376 of file Amesos_Mumps.h.

int Amesos_Mumps::MtxConvTime_
protected

Quick access pointers to the internal timers.

Definition at line 379 of file Amesos_Mumps.h.

int Amesos_Mumps::MtxRedistTime_
protected

Definition at line 379 of file Amesos_Mumps.h.

int Amesos_Mumps::VecRedistTime_
protected

Definition at line 379 of file Amesos_Mumps.h.

int Amesos_Mumps::SymFactTime_
protected

Definition at line 380 of file Amesos_Mumps.h.

int Amesos_Mumps::NumFactTime_
protected

Definition at line 380 of file Amesos_Mumps.h.

int Amesos_Mumps::SolveTime_
protected

Definition at line 380 of file Amesos_Mumps.h.

double* Amesos_Mumps::RowSca_
protected

Row and column scaling.

Definition at line 383 of file Amesos_Mumps.h.

double * Amesos_Mumps::ColSca_
protected

Definition at line 383 of file Amesos_Mumps.h.

int* Amesos_Mumps::PermIn_
protected

PermIn for MUMPS.

Definition at line 386 of file Amesos_Mumps.h.

int Amesos_Mumps::NumSchurComplementRows_
protected

Number of rows in the Schur complement (if required)

Definition at line 389 of file Amesos_Mumps.h.

int* Amesos_Mumps::SchurComplementRows_
protected

Rows for the Schur complement (if required)

Definition at line 391 of file Amesos_Mumps.h.

RCP<Epetra_CrsMatrix> Amesos_Mumps::CrsSchurComplement_
protected

Pointer to the Schur complement, as CrsMatrix.

Definition at line 394 of file Amesos_Mumps.h.

RCP<Epetra_SerialDenseMatrix> Amesos_Mumps::DenseSchurComplement_
protected

Pointer to the Schur complement,as DenseMatrix.

Definition at line 396 of file Amesos_Mumps.h.

const Epetra_LinearProblem* Amesos_Mumps::Problem_
protected

Pointer to the linear problem to be solved.

Definition at line 404 of file Amesos_Mumps.h.

RCP<Epetra_Map> Amesos_Mumps::RedistrMap_
protected

Redistributed matrix.

Definition at line 407 of file Amesos_Mumps.h.

RCP<Epetra_Import> Amesos_Mumps::RedistrImporter_
protected

Redistributed importer (from Matrix().RowMatrixRowMap() to RedistrMatrix().RowMatrixRowMap()).

Definition at line 409 of file Amesos_Mumps.h.

RCP<Epetra_CrsMatrix> Amesos_Mumps::RedistrMatrix_
protected

Redistributed matrix (only if MaxProcs_ > 1).

Definition at line 411 of file Amesos_Mumps.h.

RCP<Epetra_Map> Amesos_Mumps::SerialMap_
protected

Map with all elements on process 0 (for solution and rhs).

Definition at line 413 of file Amesos_Mumps.h.

RCP<Epetra_Import> Amesos_Mumps::SerialImporter_
protected

Importer from Matrix.OperatorDomainMap() to SerialMap_.

Definition at line 415 of file Amesos_Mumps.h.

DMUMPS_STRUC_C Amesos_Mumps::MDS
protected

Definition at line 417 of file Amesos_Mumps.h.

std::map<int, int> Amesos_Mumps::ICNTL
protected

Definition at line 419 of file Amesos_Mumps.h.

std::map<int, double> Amesos_Mumps::CNTL
protected

Definition at line 420 of file Amesos_Mumps.h.


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