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

#include <EpetraExt_HypreIJMatrix.h>

Inheritance diagram for EpetraExt_HypreIJMatrix:
Inheritance graph
[legend]

Protected Member Functions

int InitializeDefaults ()
 Set global variables to default values. More...
 
int Hypre_BoomerAMGCreate (MPI_Comm comm, HYPRE_Solver *solver)
 AMG Create passing function. More...
 
int Hypre_ParaSailsCreate (MPI_Comm comm, HYPRE_Solver *solver)
 ParaSails Create passing function. More...
 
int Hypre_EuclidCreate (MPI_Comm comm, HYPRE_Solver *solver)
 Euclid Create passing function. More...
 
int Hypre_AMSCreate (MPI_Comm comm, HYPRE_Solver *solver)
 AMS Create passing function. More...
 
int Hypre_ParCSRHybridCreate (MPI_Comm comm, HYPRE_Solver *solver)
 Hybrid Create passing function. More...
 
int Hypre_ParCSRPCGCreate (MPI_Comm comm, HYPRE_Solver *solver)
 PCG Create passing function. More...
 
int Hypre_ParCSRGMRESCreate (MPI_Comm comm, HYPRE_Solver *solver)
 GMRES Create passing function. More...
 
int Hypre_ParCSRFlexGMRESCreate (MPI_Comm comm, HYPRE_Solver *solver)
 FlexGMRES Create passing function. More...
 
int Hypre_ParCSRLGMRESCreate (MPI_Comm comm, HYPRE_Solver *solver)
 LGMRES Create passing function. More...
 
int Hypre_ParCSRBiCGSTABCreate (MPI_Comm comm, HYPRE_Solver *solver)
 BiCGSTAB Create passing function. More...
 
int CreateSolver ()
 Create the solver with selected type. More...
 
int CreatePrecond ()
 Create the preconditioner with selected type. More...
 
int SetupSolver () const
 
int SetupPrecond () const
 
- Protected Member Functions inherited from Epetra_BasicRowMatrix
void Setup ()
 
void UpdateImportVector (int NumVectors) const
 
void UpdateExportVector (int NumVectors) const
 
void SetImportExport ()
 
std::string toString (const int &x) const
 
std::string toString (const long long &x) const
 
std::string toString (const double &x) const
 
virtual void ComputeStructureConstants () const
 
virtual void ComputeNumericConstants () const
 

Protected Attributes

HYPRE_IJMatrix Matrix_
 
HYPRE_ParCSRMatrix ParMatrix_
 
HYPRE_IJVector X_hypre
 
HYPRE_IJVector Y_hypre
 
HYPRE_ParVector par_x
 
HYPRE_ParVector par_y
 
hypre_ParVector * x_vec
 
hypre_ParVector * y_vec
 
hypre_Vector * x_local
 
hypre_Vector * y_local
 
HYPRE_Solver Solver_
 
HYPRE_Solver Preconditioner_
 
int(EpetraExt_HypreIJMatrix::* SolverCreatePtr_ )(MPI_Comm, HYPRE_Solver *)
 
int(* SolverDestroyPtr_ )(HYPRE_Solver)
 
int(* SolverSetupPtr_ )(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
 
int(* SolverSolvePtr_ )(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
 
int(* SolverPrecondPtr_ )(HYPRE_Solver, HYPRE_PtrToParSolverFcn, HYPRE_PtrToParSolverFcn, HYPRE_Solver)
 
int(EpetraExt_HypreIJMatrix::* PrecondCreatePtr_ )(MPI_Comm, HYPRE_Solver *)
 
int(* PrecondDestroyPtr_ )(HYPRE_Solver)
 
int(* PrecondSetupPtr_ )(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
 
int(* PrecondSolvePtr_ )(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
 
int NumMyRows_
 Number of rows on local processor. More...
 
int NumGlobalRows_
 Number of rows across all processors. More...
 
int NumGlobalCols_
 Number of columns across all processors. More...
 
int MyRowStart_
 First row on local processor. More...
 
int MyRowEnd_
 Last row on local processor. More...
 
int MatType_
 Hypre matrix type (parCSR). More...
 
bool TransposeSolve_
 Do a transpose solve, only BoomerAMG. More...
 
Hypre_Chooser SolveOrPrec_
 Choose to solve or apply preconditioner. More...
 
bool * IsSolverSetup_
 Flag to know if solver needs to be destoyed. More...
 
bool * IsPrecondSetup_
 Flag to know if preconditioner needs to be destroyed. More...
 
- Protected Attributes inherited from Epetra_BasicRowMatrix
Epetra_CommComm_
 
Epetra_Map OperatorDomainMap_
 
Epetra_Map OperatorRangeMap_
 
Epetra_Map RowMatrixRowMap_
 
Epetra_Map RowMatrixColMap_
 
int NumMyNonzeros_
 
long long NumGlobalNonzeros_
 
int MaxNumEntries_
 
double NormInf_
 
double NormOne_
 
int NumMyRows_
 
int NumMyCols_
 
bool UseTranspose_
 
bool HasNormInf_
 
bool LowerTriangular_
 
bool UpperTriangular_
 
bool HaveStructureConstants_
 
bool HaveNumericConstants_
 
bool HaveMaps_
 
Epetra_MultiVectorImportVector_
 
Epetra_MultiVectorExportVector_
 
Epetra_ImportImporter_
 
Epetra_ExportExporter_
 
Epetra_FlopsFlopCounter_
 

Constructors/Destructor

 EpetraExt_HypreIJMatrix (HYPRE_IJMatrix matrix)
 Epetra_HypreIJMatrix constructor. More...
 
virtual ~EpetraExt_HypreIJMatrix ()
 EpetraExt_HypreIJMatrix Destructor. More...
 

Extraction methods

int ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const
 Returns a copy of the specified local row in user-provided arrays. More...
 
int ExtractMyEntryView (int CurEntry, double *&Value, int &RowIndex, int &ColIndex)
 Returns a reference to the ith entry in the matrix, along with its row and column index. More...
 
int ExtractMyEntryView (int CurEntry, const double *&Value, int &RowIndex, int &ColIndex) const
 Returns a const reference to the ith entry in the matrix, along with its row and column index. More...
 

Solver Setup Methods

int SetParameter (Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, int), int parameter)
 Set a parameter that takes a single int. More...
 
int SetParameter (Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, double), double parameter)
 Set a parameter that takes a single double. More...
 
int SetParameter (Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, double, int), double parameter1, int parameter2)
 Set a parameter that takes a double then an int. More...
 
int SetParameter (Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, int, int), int parameter1, int parameter2)
 Set a parameter that takes two int parameters. More...
 
int SetParameter (Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, double *), double *parameter)
 Set a parameter that takes a double*. More...
 
int SetParameter (Hypre_Chooser chooser, int(*pt2Func)(HYPRE_Solver, int *), int *parameter)
 Set a parameter that takes an int*. More...
 
int SetParameter (Hypre_Chooser chooser, Hypre_Solver Solver, bool transpose=false)
 Sets the solver that is used by the Solve() and ApplyInverse() methods. Until this is called, the default solver is PCG. More...
 
int SetParameter (bool UsePreconditioner)
 Sets the solver to use the selected preconditioner. More...
 
int SetParameter (Hypre_Chooser chooser)
 Choose to solve the problem or apply the preconditioner. More...
 

Computational methods

int Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a EpetraExt_HypreIJMatrix multiplied by a Epetra_MultiVector X in Y. More...
 
int Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a EpetraExt_HypreIJMatrix solving a Epetra_MultiVector X in Y. More...
 
int LeftScale (const Epetra_Vector &X)
 Scales the EpetraExt_HypreIJMatrix on the left with a Epetra_Vector x. More...
 
int RightScale (const Epetra_Vector &X)
 Scales the EpetraExt_HypreIJMatrix on the right with a Epetra_Vector x. More...
 

Additional methods required to support the Epetra_Operator interface

int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y. More...
 
int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y. More...
 
virtual bool UseTranspose () const
 Returns the current UseTranspose setting. More...
 

Additional methods required to implement RowMatrix interface

int NumMyRowEntries (int MyRow, int &NumEntries) const
 Return the current number of values stored for the specified local row. More...
 
HYPRE_IJMatrix & GetMatrix ()
 Return a reference to the Hypre matrix. More...
 

Additional Inherited Members

- Public Member Functions inherited from Epetra_BasicRowMatrix
 Epetra_BasicRowMatrix (const Epetra_Comm &Comm)
 
virtual ~Epetra_BasicRowMatrix ()
 
void SetMaps (const Epetra_Map &RowMap, const Epetra_Map &ColMap)
 
void SetMaps (const Epetra_Map &RowMap, const Epetra_Map &ColMap, const Epetra_Map &DomainMap, const Epetra_Map &RangeMap)
 
virtual int ExtractDiagonalCopy (Epetra_Vector &Diagonal) const
 
virtual int InvRowSums (Epetra_Vector &x) const
 
virtual int InvColSums (Epetra_Vector &x) const
 
virtual bool Filled () const
 
bool LowerTriangular () const
 
virtual bool UpperTriangular () const
 
virtual double NormInf () const
 
virtual double NormOne () const
 
virtual int NumGlobalNonzeros () const
 
virtual long long NumGlobalNonzeros64 () const
 
virtual int NumGlobalRows () const
 
virtual long long NumGlobalRows64 () const
 
virtual int NumGlobalCols () const
 
virtual long long NumGlobalCols64 () const
 
virtual int NumGlobalDiagonals () const
 
virtual long long NumGlobalDiagonals64 () const
 
virtual int NumMyNonzeros () const
 
virtual int NumMyRows () const
 
virtual int NumMyCols () const
 
virtual int NumMyDiagonals () const
 
virtual int MaxNumEntries () const
 
virtual const Epetra_MapOperatorDomainMap () const
 
virtual const Epetra_MapOperatorRangeMap () const
 
virtual const Epetra_BlockMapMap () const
 
virtual const Epetra_MapRowMatrixRowMap () const
 
virtual const Epetra_MapRowMatrixColMap () const
 
virtual const Epetra_ImportRowMatrixImporter () const
 
virtual const Epetra_CommComm () const
 
virtual void Print (std::ostream &os) const
 
virtual int SetUseTranspose (bool use_transpose)
 
virtual const char * Label () const
 
bool HasNormInf () const
 
virtual const Epetra_ImportImporter () const
 
virtual const Epetra_ExportExporter () const
 
Epetra_CompObjectoperator= (const Epetra_CompObject &src)
 
 Epetra_CompObject ()
 
 Epetra_CompObject (const Epetra_CompObject &Source)
 
virtual ~Epetra_CompObject ()
 
void SetFlopCounter (const Epetra_Flops &FlopCounter_in)
 
void SetFlopCounter (const Epetra_CompObject &CompObject)
 
void UnsetFlopCounter ()
 
Epetra_FlopsGetFlopCounter () const
 
void ResetFlops () const
 
double Flops () const
 
void UpdateFlops (int Flops_in) const
 
void UpdateFlops (long int Flops_in) const
 
void UpdateFlops (long long Flops_in) const
 
void UpdateFlops (double Flops_in) const
 
void UpdateFlops (float Flops_in) const
 
 Epetra_Object (int TracebackModeIn=-1, bool set_label=true)
 
 Epetra_Object (const char *const Label, int TracebackModeIn=-1)
 
 Epetra_Object (const Epetra_Object &Object)
 
virtual ~Epetra_Object ()
 
virtual void SetLabel (const char *const Label)
 
virtual int ReportError (const std::string Message, int ErrorCode) const
 
virtual ~Epetra_RowMatrix ()
 
virtual ~Epetra_Operator ()
 
virtual ~Epetra_SrcDistObject ()
 
- Static Public Member Functions inherited from Epetra_BasicRowMatrix
static void SetTracebackMode (int TracebackModeValue)
 
static int GetTracebackMode ()
 
static std::ostream & GetTracebackStream ()
 
- Static Public Attributes inherited from Epetra_BasicRowMatrix
static int TracebackMode
 

Detailed Description

Definition at line 99 of file EpetraExt_HypreIJMatrix.h.

Constructor & Destructor Documentation

EpetraExt_HypreIJMatrix::EpetraExt_HypreIJMatrix ( HYPRE_IJMatrix  matrix)

Epetra_HypreIJMatrix constructor.

Creates a Epetra_HypreIJMatrix object by encapsulating an existing Hypre matrix.

Parameters
matrix(In) - A completely constructed Hypre IJ matrix.
virtual EpetraExt_HypreIJMatrix::~EpetraExt_HypreIJMatrix ( )
virtual

Member Function Documentation

int EpetraExt_HypreIJMatrix::ExtractMyRowCopy ( int  MyRow,
int  Length,
int &  NumEntries,
double *  Values,
int *  Indices 
) const
virtual

Returns a copy of the specified local row in user-provided arrays.

Parameters
MyRow(In) - Local row to extract.
Length(In) - Length of Values and Indices.
NumEntries(Out) - Number of nonzero entries extracted.
Values(Out) - Extracted values for this row.
Indices(Out) - Extracted local column indices for the corresponding values.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_BasicRowMatrix.

int EpetraExt_HypreIJMatrix::ExtractMyEntryView ( int  CurEntry,
double *&  Value,
int &  RowIndex,
int &  ColIndex 
)
virtual

Returns a reference to the ith entry in the matrix, along with its row and column index.

Parameters
CurEntry(In) - Index of local entry (from 0 to NumMyNonzeros()-1) to extract.
Value(Out) - Extracted reference to current values.
RowIndex(Out) - Row index for current entry.
ColIndex(Out) - Column index for current entry.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_BasicRowMatrix.

int EpetraExt_HypreIJMatrix::ExtractMyEntryView ( int  CurEntry,
const double *&  Value,
int &  RowIndex,
int &  ColIndex 
) const
virtual

Returns a const reference to the ith entry in the matrix, along with its row and column index.

Parameters
CurEntry(In) - Index of local entry (from 0 to NumMyNonzeros()-1) to extract.
Value(Out) - Extracted reference to current values.
RowIndex(Out) - Row index for current entry.
ColIndex(Out) - Column index for current entry.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_BasicRowMatrix.

int EpetraExt_HypreIJMatrix::SetParameter ( Hypre_Chooser  chooser,
int(*)(HYPRE_Solver, int)  pt2Func,
int  parameter 
)

Set a parameter that takes a single int.

Parameters
chooser(In) -A Hypre_Chooser enumerated type set to Solver or Preconditioner, whatever the parameter is setting for.
*pt2Func(In) -The function that sets the parameter. It must set parameters for the type of solver or preconditioner that was created. An example is if the solver is BoomerAMG, the function to set maximum iterations would be &HYPRE_BoomerAMGSetMaxIter
parameter(In) -The integer parameter being set.
Returns
Integer error code, set to 0 if successful.
int EpetraExt_HypreIJMatrix::SetParameter ( Hypre_Chooser  chooser,
int(*)(HYPRE_Solver, double)  pt2Func,
double  parameter 
)

Set a parameter that takes a single double.

Parameters
chooser(In) -A Hypre_Chooser enumerated type set to Solver or Preconditioner, whatever the parameter is setting for.
*pt2Func(In) -The function that sets the parameter. It must set parameters for the type of solver or preconditioner that was created. An example is if the solver is BoomerAMG, the function to set tolerance would be &HYPRE_BoomerAMGSetTol
parameter(In) -The double parameter being set.
Returns
Integer error code, set to 0 if successful.
int EpetraExt_HypreIJMatrix::SetParameter ( Hypre_Chooser  chooser,
int(*)(HYPRE_Solver, double, int)  pt2Func,
double  parameter1,
int  parameter2 
)

Set a parameter that takes a double then an int.

Parameters
chooser(In) -A Hypre_Chooser enumerated type set to Solver or Preconditioner, whatever the parameter is setting for.
*pt2Func(In) -The function that sets the parameter. It must set parameters for the type of solver or preconditioner that was created. An example is if the solver is BoomerAMG, the function to set relaxation weight for a given level would be &HYPRE_BoomerAMGSetLevelRelaxWt
parameter1(In) -The double parameter being set.
parameter2(In) - The integer parameter being set.
Returns
Integer error code, set to 0 if successful.
int EpetraExt_HypreIJMatrix::SetParameter ( Hypre_Chooser  chooser,
int(*)(HYPRE_Solver, int, int)  pt2Func,
int  parameter1,
int  parameter2 
)

Set a parameter that takes two int parameters.

Parameters
chooser(In) -A Hypre_Chooser enumerated type set to Solver or Preconditioner, whatever the parameter is setting for.
*pt2Func(In) -The function that sets the parameter. It must set parameters for the type of solver or preconditioner that was created. An example is if the solver is BoomerAMG, the function to set relaxation type for a given level would be &HYPRE_BoomerAMGSetCycleRelaxType
parameter1(In) -The first integer parameter being set.
parameter2(In) - The second integer parameter being set.
Returns
Integer error code, set to 0 if successful.
int EpetraExt_HypreIJMatrix::SetParameter ( Hypre_Chooser  chooser,
int(*)(HYPRE_Solver, double *)  pt2Func,
double *  parameter 
)

Set a parameter that takes a double*.

Parameters
chooser(In) -A Hypre_Chooser enumerated type set to Solver or Preconditioner, whatever the parameter is setting for.
*pt2Func(In) -The function that sets the parameter. It must set parameters for the type of solver or preconditioner that was created. An example is if the solver is BoomerAMG, the function to set relaxation weight would be &HYPRE_BoomerAMGSetRelaxWeight
parameter(In) -The double* parameter being set.
Returns
Integer error code, set to 0 if successful.
int EpetraExt_HypreIJMatrix::SetParameter ( Hypre_Chooser  chooser,
int(*)(HYPRE_Solver, int *)  pt2Func,
int *  parameter 
)

Set a parameter that takes an int*.

Parameters
chooser(In) -A Hypre_Chooser enumerated type set to Solver or Preconditioner, whatever the parameter is setting for.
*pt2Func(In) -The function that sets the parameter. It must set parameters for the type of solver or preconditioner that was created. An example is if the solver is BoomerAMG, the function to set grid relax type would be &HYPRE_BoomerAMGSetGridRelaxType
parameter(In) -The int* parameter being set.
Returns
Integer error code, set to 0 if successful.
int EpetraExt_HypreIJMatrix::SetParameter ( Hypre_Chooser  chooser,
Hypre_Solver  Solver,
bool  transpose = false 
)

Sets the solver that is used by the Solve() and ApplyInverse() methods. Until this is called, the default solver is PCG.

Parameters
chooser(In) - A Hypre_Chooser enumerated type. If Solver, then we are selecting which solver, if Preconditioner, we are choosing which preconditioner to use.
Solver(In) -A Hypre_Solver enumerated type to select the solver or preconditioner. Options for solver are: BoomerAMG, AMS, Hybrid, PCG, GMRES, FlexGMRES, LGMRES, and BiCGSTAB. See Hypre Ref Manual for more info on the solvers. Options for Preconditioner are: BoomerAMG, ParaSails, Euclid, and AMS.
transpose(In) -Optional argument that selects to use a transpose solve. Currently BoomerAMG is the only solver available for transpose solve. It must be the argument for Solver if transpose is true.
Returns
Integer error code, set to 0 if successful.
int EpetraExt_HypreIJMatrix::SetParameter ( bool  UsePreconditioner)

Sets the solver to use the selected preconditioner.

Parameters
UsePreconditioner(In) -A boolean, true use preconditioner, false do not use the supplied preconditioner with the solver. The solver and preconditioner must have been selected and the solver must be one of the following solvers: Hybrid, PCG, GMRES, FlexGMRES, LGMRES, BiCGSTAB.
Returns
Integer error code, set to 0 if successful.
int EpetraExt_HypreIJMatrix::SetParameter ( Hypre_Chooser  chooser)

Choose to solve the problem or apply the preconditioner.

Parameters
chooser(In) -A Hypre_Chooser enumerated type, either Solver or Preconditioner. The chosen type must have been selected before this method is called.
Returns
Integer error code, set to 0 if successful.
int EpetraExt_HypreIJMatrix::Multiply ( bool  TransA,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Returns the result of a EpetraExt_HypreIJMatrix multiplied by a Epetra_MultiVector X in Y.

Parameters
TransA(In) -If true, multiply by the transpose of matrix, otherwise just use matrix.
X(In) - A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
Y(Out) -A Epetra_MultiVector of dimension NumVectorscontaining result.
Returns
Integer error code, set to 0 if successful.

Reimplemented from Epetra_BasicRowMatrix.

int EpetraExt_HypreIJMatrix::Solve ( bool  Upper,
bool  Trans,
bool  UnitDiagonal,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

Returns the result of a EpetraExt_HypreIJMatrix solving a Epetra_MultiVector X in Y.

Parameters
Upper(In) -If true, solve Ux = y, otherwise solve Lx = y.
Trans(In) -If true, solve transpose problem.
UnitDiagonal(In) -If true, assume diagonal is unit (whether it's stored or not).
X(In) - A Epetra_MultiVector of dimension NumVectors to solve for.
Y(Out) -A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Reimplemented from Epetra_BasicRowMatrix.

int EpetraExt_HypreIJMatrix::LeftScale ( const Epetra_Vector X)
virtual

Scales the EpetraExt_HypreIJMatrix on the left with a Epetra_Vector x.

The this matrix will be scaled such that A(i,j) = x(i)*A(i,j) where i denotes the row number of A and j denotes the column number of A.

Parameters
X(In) -A Epetra_Vector to solve for.
Returns
Integer error code, set to 0 if successful.

Reimplemented from Epetra_BasicRowMatrix.

int EpetraExt_HypreIJMatrix::RightScale ( const Epetra_Vector X)
virtual

Scales the EpetraExt_HypreIJMatrix on the right with a Epetra_Vector x.

The this matrix will be scaled such that A(i,j) = x(j)*A(i,j) where i denotes the global row number of A and j denotes the global column number of A.

Parameters
X(In) -The Epetra_Vector used for scaling this.
Returns
Integer error code, set to 0 if successful.

Reimplemented from Epetra_BasicRowMatrix.

int EpetraExt_HypreIJMatrix::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.

Parameters
X(In) - A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
Y(Out) -A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Reimplemented from Epetra_BasicRowMatrix.

Definition at line 316 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.

In this implementation, we use several existing attributes to determine how virtual
method ApplyInverse() should call the concrete method Solve().  We pass in the UpperTriangular(), 

the EpetraExt_HypreIJMatrix::UseTranspose(), and NoDiagonal() methods. The most notable warning is that if a matrix has no diagonal values we assume that there is an implicit unit diagonal that should be accounted for when doing a triangular solve.

Parameters
X(In) - A Epetra_MultiVector of dimension NumVectors to solve for.
Y(Out) -A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Reimplemented from Epetra_BasicRowMatrix.

Definition at line 331 of file EpetraExt_HypreIJMatrix.h.

virtual bool EpetraExt_HypreIJMatrix::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Reimplemented from Epetra_BasicRowMatrix.

Definition at line 335 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::NumMyRowEntries ( int  MyRow,
int &  NumEntries 
) const
virtual

Return the current number of values stored for the specified local row.

Similar to NumMyEntries() except NumEntries is returned as an argument
and error checking is done on the input value MyRow.
Parameters
MyRow(In) - Local row.
NumEntries(Out) - Number of nonzero values.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_BasicRowMatrix.

HYPRE_IJMatrix& EpetraExt_HypreIJMatrix::GetMatrix ( )
inline

Return a reference to the Hypre matrix.

Definition at line 352 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::InitializeDefaults ( )
protected

Set global variables to default values.

int EpetraExt_HypreIJMatrix::Hypre_BoomerAMGCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

AMG Create passing function.

Definition at line 363 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::Hypre_ParaSailsCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

ParaSails Create passing function.

Definition at line 367 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::Hypre_EuclidCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

Euclid Create passing function.

Definition at line 371 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::Hypre_AMSCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

AMS Create passing function.

Definition at line 375 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::Hypre_ParCSRHybridCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

Hybrid Create passing function.

Definition at line 379 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::Hypre_ParCSRPCGCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

PCG Create passing function.

Definition at line 383 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::Hypre_ParCSRGMRESCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

GMRES Create passing function.

Definition at line 387 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::Hypre_ParCSRFlexGMRESCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

FlexGMRES Create passing function.

Definition at line 391 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::Hypre_ParCSRLGMRESCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

LGMRES Create passing function.

Definition at line 395 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::Hypre_ParCSRBiCGSTABCreate ( MPI_Comm  comm,
HYPRE_Solver *  solver 
)
inlineprotected

BiCGSTAB Create passing function.

Definition at line 399 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::CreateSolver ( )
protected

Create the solver with selected type.

int EpetraExt_HypreIJMatrix::CreatePrecond ( )
protected

Create the preconditioner with selected type.

int EpetraExt_HypreIJMatrix::SetupSolver ( ) const
protected
int EpetraExt_HypreIJMatrix::SetupPrecond ( ) const
protected

Member Data Documentation

HYPRE_IJMatrix EpetraExt_HypreIJMatrix::Matrix_
mutableprotected

Definition at line 413 of file EpetraExt_HypreIJMatrix.h.

HYPRE_ParCSRMatrix EpetraExt_HypreIJMatrix::ParMatrix_
mutableprotected

Definition at line 414 of file EpetraExt_HypreIJMatrix.h.

HYPRE_IJVector EpetraExt_HypreIJMatrix::X_hypre
mutableprotected

Definition at line 415 of file EpetraExt_HypreIJMatrix.h.

HYPRE_IJVector EpetraExt_HypreIJMatrix::Y_hypre
mutableprotected

Definition at line 416 of file EpetraExt_HypreIJMatrix.h.

HYPRE_ParVector EpetraExt_HypreIJMatrix::par_x
mutableprotected

Definition at line 417 of file EpetraExt_HypreIJMatrix.h.

HYPRE_ParVector EpetraExt_HypreIJMatrix::par_y
mutableprotected

Definition at line 418 of file EpetraExt_HypreIJMatrix.h.

hypre_ParVector* EpetraExt_HypreIJMatrix::x_vec
mutableprotected

Definition at line 419 of file EpetraExt_HypreIJMatrix.h.

hypre_ParVector* EpetraExt_HypreIJMatrix::y_vec
mutableprotected

Definition at line 420 of file EpetraExt_HypreIJMatrix.h.

hypre_Vector* EpetraExt_HypreIJMatrix::x_local
mutableprotected

Definition at line 421 of file EpetraExt_HypreIJMatrix.h.

hypre_Vector* EpetraExt_HypreIJMatrix::y_local
mutableprotected

Definition at line 422 of file EpetraExt_HypreIJMatrix.h.

HYPRE_Solver EpetraExt_HypreIJMatrix::Solver_
mutableprotected

Definition at line 423 of file EpetraExt_HypreIJMatrix.h.

HYPRE_Solver EpetraExt_HypreIJMatrix::Preconditioner_
mutableprotected

Definition at line 424 of file EpetraExt_HypreIJMatrix.h.

int(EpetraExt_HypreIJMatrix::* EpetraExt_HypreIJMatrix::SolverCreatePtr_)(MPI_Comm, HYPRE_Solver *)
protected

Definition at line 426 of file EpetraExt_HypreIJMatrix.h.

int(* EpetraExt_HypreIJMatrix::SolverDestroyPtr_)(HYPRE_Solver)
protected

Definition at line 427 of file EpetraExt_HypreIJMatrix.h.

int(* EpetraExt_HypreIJMatrix::SolverSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
protected

Definition at line 428 of file EpetraExt_HypreIJMatrix.h.

int(* EpetraExt_HypreIJMatrix::SolverSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
protected

Definition at line 429 of file EpetraExt_HypreIJMatrix.h.

int(* EpetraExt_HypreIJMatrix::SolverPrecondPtr_)(HYPRE_Solver, HYPRE_PtrToParSolverFcn, HYPRE_PtrToParSolverFcn, HYPRE_Solver)
protected

Definition at line 430 of file EpetraExt_HypreIJMatrix.h.

int(EpetraExt_HypreIJMatrix::* EpetraExt_HypreIJMatrix::PrecondCreatePtr_)(MPI_Comm, HYPRE_Solver *)
protected

Definition at line 431 of file EpetraExt_HypreIJMatrix.h.

int(* EpetraExt_HypreIJMatrix::PrecondDestroyPtr_)(HYPRE_Solver)
protected

Definition at line 432 of file EpetraExt_HypreIJMatrix.h.

int(* EpetraExt_HypreIJMatrix::PrecondSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
protected

Definition at line 433 of file EpetraExt_HypreIJMatrix.h.

int(* EpetraExt_HypreIJMatrix::PrecondSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector)
protected

Definition at line 434 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::NumMyRows_
protected

Number of rows on local processor.

Definition at line 437 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::NumGlobalRows_
protected

Number of rows across all processors.

Definition at line 439 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::NumGlobalCols_
protected

Number of columns across all processors.

Definition at line 441 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::MyRowStart_
protected

First row on local processor.

Definition at line 443 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::MyRowEnd_
protected

Last row on local processor.

Definition at line 445 of file EpetraExt_HypreIJMatrix.h.

int EpetraExt_HypreIJMatrix::MatType_
mutableprotected

Hypre matrix type (parCSR).

Definition at line 447 of file EpetraExt_HypreIJMatrix.h.

bool EpetraExt_HypreIJMatrix::TransposeSolve_
protected

Do a transpose solve, only BoomerAMG.

Definition at line 449 of file EpetraExt_HypreIJMatrix.h.

Hypre_Chooser EpetraExt_HypreIJMatrix::SolveOrPrec_
protected

Choose to solve or apply preconditioner.

Definition at line 451 of file EpetraExt_HypreIJMatrix.h.

bool* EpetraExt_HypreIJMatrix::IsSolverSetup_
protected

Flag to know if solver needs to be destoyed.

Definition at line 453 of file EpetraExt_HypreIJMatrix.h.

bool* EpetraExt_HypreIJMatrix::IsPrecondSetup_
protected

Flag to know if preconditioner needs to be destroyed.

Definition at line 455 of file EpetraExt_HypreIJMatrix.h.


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