#include <EpetraExt_HypreIJMatrix.h>
|
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 |
|
virtual void | ComputeStructureConstants () const |
|
virtual void | ComputeNumericConstants () const |
|
|
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...
|
|
|
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...
|
|
|
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...
|
|
|
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...
|
|
Definition at line 99 of file EpetraExt_HypreIJMatrix.h.
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 |
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.
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.
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
-
- 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
-
- Returns
- Integer error code, set to 0 if successful.
Reimplemented from Epetra_BasicRowMatrix.
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
-
- 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 |
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 |
int EpetraExt_HypreIJMatrix::InitializeDefaults |
( |
| ) |
|
|
protected |
Set global variables to default values.
int EpetraExt_HypreIJMatrix::Hypre_BoomerAMGCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
int EpetraExt_HypreIJMatrix::Hypre_ParaSailsCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
int EpetraExt_HypreIJMatrix::Hypre_EuclidCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
int EpetraExt_HypreIJMatrix::Hypre_AMSCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
int EpetraExt_HypreIJMatrix::Hypre_ParCSRHybridCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
int EpetraExt_HypreIJMatrix::Hypre_ParCSRPCGCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
int EpetraExt_HypreIJMatrix::Hypre_ParCSRGMRESCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
int EpetraExt_HypreIJMatrix::Hypre_ParCSRFlexGMRESCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
int EpetraExt_HypreIJMatrix::Hypre_ParCSRLGMRESCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
int EpetraExt_HypreIJMatrix::Hypre_ParCSRBiCGSTABCreate |
( |
MPI_Comm |
comm, |
|
|
HYPRE_Solver * |
solver |
|
) |
| |
|
inlineprotected |
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 |
HYPRE_IJMatrix EpetraExt_HypreIJMatrix::Matrix_ |
|
mutableprotected |
HYPRE_ParCSRMatrix EpetraExt_HypreIJMatrix::ParMatrix_ |
|
mutableprotected |
HYPRE_IJVector EpetraExt_HypreIJMatrix::X_hypre |
|
mutableprotected |
HYPRE_IJVector EpetraExt_HypreIJMatrix::Y_hypre |
|
mutableprotected |
HYPRE_ParVector EpetraExt_HypreIJMatrix::par_x |
|
mutableprotected |
HYPRE_ParVector EpetraExt_HypreIJMatrix::par_y |
|
mutableprotected |
hypre_ParVector* EpetraExt_HypreIJMatrix::x_vec |
|
mutableprotected |
hypre_ParVector* EpetraExt_HypreIJMatrix::y_vec |
|
mutableprotected |
hypre_Vector* EpetraExt_HypreIJMatrix::x_local |
|
mutableprotected |
hypre_Vector* EpetraExt_HypreIJMatrix::y_local |
|
mutableprotected |
HYPRE_Solver EpetraExt_HypreIJMatrix::Solver_ |
|
mutableprotected |
HYPRE_Solver EpetraExt_HypreIJMatrix::Preconditioner_ |
|
mutableprotected |
int(EpetraExt_HypreIJMatrix::* EpetraExt_HypreIJMatrix::SolverCreatePtr_)(MPI_Comm, HYPRE_Solver *) |
|
protected |
int(* EpetraExt_HypreIJMatrix::SolverDestroyPtr_)(HYPRE_Solver) |
|
protected |
int(* EpetraExt_HypreIJMatrix::SolverSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector) |
|
protected |
int(* EpetraExt_HypreIJMatrix::SolverSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector) |
|
protected |
int(* EpetraExt_HypreIJMatrix::SolverPrecondPtr_)(HYPRE_Solver, HYPRE_PtrToParSolverFcn, HYPRE_PtrToParSolverFcn, HYPRE_Solver) |
|
protected |
int(EpetraExt_HypreIJMatrix::* EpetraExt_HypreIJMatrix::PrecondCreatePtr_)(MPI_Comm, HYPRE_Solver *) |
|
protected |
int(* EpetraExt_HypreIJMatrix::PrecondDestroyPtr_)(HYPRE_Solver) |
|
protected |
int(* EpetraExt_HypreIJMatrix::PrecondSetupPtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector) |
|
protected |
int(* EpetraExt_HypreIJMatrix::PrecondSolvePtr_)(HYPRE_Solver, HYPRE_ParCSRMatrix, HYPRE_ParVector, HYPRE_ParVector) |
|
protected |
int EpetraExt_HypreIJMatrix::NumMyRows_ |
|
protected |
int EpetraExt_HypreIJMatrix::NumGlobalRows_ |
|
protected |
int EpetraExt_HypreIJMatrix::NumGlobalCols_ |
|
protected |
int EpetraExt_HypreIJMatrix::MyRowStart_ |
|
protected |
int EpetraExt_HypreIJMatrix::MyRowEnd_ |
|
protected |
int EpetraExt_HypreIJMatrix::MatType_ |
|
mutableprotected |
bool EpetraExt_HypreIJMatrix::TransposeSolve_ |
|
protected |
bool* EpetraExt_HypreIJMatrix::IsSolverSetup_ |
|
protected |
bool* EpetraExt_HypreIJMatrix::IsPrecondSetup_ |
|
protected |
The documentation for this class was generated from the following file: