Pliris: An Obect-Oriented Interface to a Dense LU Solver. More...
#include <Pliris.h>
Public Member Functions | |
Pliris (Epetra_Vector *A, Epetra_MultiVector *X, Epetra_MultiVector *B) | |
Pliris () | |
Pliris Default constructor. | |
int | SetLHS (Epetra_MultiVector *X) |
Pliris LHS Set. More... | |
int | SetRHS (Epetra_MultiVector *B) |
Pliris RHS Set. More... | |
int | SetMatrix (Epetra_Vector *A) |
Pliris Matrix Set. More... | |
int | SetMatrix (Epetra_SerialDenseVector *A) |
Pliris Matrix Set. More... | |
int | GetDistribution (int *nprocs_row, int *number_of_unknowns, int *nrhs, int *my_rows, int *my_cols, int *my_first_row, int *my_first_col, int *my_rhs, int *my_row, int *my_col) |
Pliris GetDistirbution. More... | |
int | FactorSolve (Epetra_Vector *A, int my_rows, int my_cols, int *matrix_size, int *num_procsr, int *num_rhs, double *secs) |
Pliris FactorSolve. More... | |
int | FactorSolve (Epetra_SerialDenseVector *AA, int my_rows, int my_cols, int *matrix_size, int *num_procsr, int *num_rhs, double *secs) |
Pliris FactorSolve. More... | |
int | Factor (Epetra_Vector *A, int *matrix_size, int *num_procsr, int *permute, double *secs) |
Pliris Factor. More... | |
int | Solve (int *permute, int *num_rhs) |
Pliris Solve. More... | |
virtual | ~Pliris (void) |
Pliris Default Destructor. More... | |
Protected Attributes | |
double * | x_ |
double * | a_ |
double * | b_ |
int | x_LDA_ |
int | b_LDA_ |
bool | inConstructor_ |
Epetra_MultiVector * | X_ |
Epetra_MultiVector * | B_ |
Epetra_Vector * | A_ |
Epetra_SerialDenseVector * | AA_ |
Pliris: An Obect-Oriented Interface to a Dense LU Solver.
The Pliris class : Provides the functionality to interface to a dense LU
int Pliris::Factor | ( | Epetra_Vector * | A, |
int * | matrix_size, | ||
int * | num_procsr, | ||
int * | permute, | ||
double * | secs | ||
) |
Pliris Factor.
Factors the dense matrix
A(In) | – Epetra Vector that has the matrix packed |
matrix_size(In) | – order of the dense matrix |
num_procsr(In) | – number of processors for a row |
permute(In) | – permutation matrix |
secs(Out) | – factor and solve time in seconds |
References SetMatrix().
int Pliris::FactorSolve | ( | Epetra_Vector * | A, |
int | my_rows, | ||
int | my_cols, | ||
int * | matrix_size, | ||
int * | num_procsr, | ||
int * | num_rhs, | ||
double * | secs | ||
) |
Pliris FactorSolve.
Factors and solves the dense matrix
A(InOut) | – Epetra Vector that has the matrix and rhs packed( Note: matrix is overwritten) |
my_rows(In) | – number of rows of the matrix on this processor |
my_cols(In) | – number of columns of the matrix on this processor |
matrix_size(In) | – order of the dense matrix |
num_procsr(In) | – number of processors for a row |
num_rhs(In) | – number of right hand sides |
secs(Out) | – factor and solve time in seconds |
References SetMatrix().
int Pliris::FactorSolve | ( | Epetra_SerialDenseVector * | AA, |
int | my_rows, | ||
int | my_cols, | ||
int * | matrix_size, | ||
int * | num_procsr, | ||
int * | num_rhs, | ||
double * | secs | ||
) |
Pliris FactorSolve.
Factors and solves the dense matrix
AA(In) | – Epetra Serial Dense Vector that has the matrix and rhs packed |
my_rows(In) | – number of rows of the matrix on this processor |
my_cols(In) | – number of columns of the matrix on this processor |
matrix_size(In) | – order of the dense matrix |
num_procsr(In) | – number of processors for a row |
num_rhs(In) | – number of right hand sides |
secs(Out) | – factor and solve time in seconds |
References SetMatrix().
int Pliris::GetDistribution | ( | int * | nprocs_row, |
int * | number_of_unknowns, | ||
int * | nrhs, | ||
int * | my_rows, | ||
int * | my_cols, | ||
int * | my_first_row, | ||
int * | my_first_col, | ||
int * | my_rhs, | ||
int * | my_row, | ||
int * | my_col | ||
) |
Pliris GetDistirbution.
Gives the distribution information that is required by the dense solver
nprocs_row(In) | - number of processors for a row |
number_of_unknowns(In) | - order of the dense matrix |
nrhs(In) | - number of right hand sides |
my_rows(Out) | - number of rows of the matrix on this processor |
my_cols | (Out) - number of columns of the matrix on this processor |
my_first_row(Out) | - first (global) row number on this processor (array starts at index 1) |
my_first_col | (Out) - first (global) column number on this processor (array starts at index 1) |
my_rhs(Out) | - number of Right hand sides on this processor |
my_row(Out) | - row number in processor mesh, 0 to the number of processors for a column -1 |
my_col(Out) | - column number in processor mesh, 0 to the number of processors for a row -1 |
int Pliris::SetLHS | ( | Epetra_MultiVector * | X | ) |
Pliris LHS Set.
Associates an already defined Epetra_MultiVector (or Epetra_Vector) as the initial guess and location where the solution will be return.
int Pliris::SetMatrix | ( | Epetra_Vector * | A | ) |
Pliris Matrix Set.
Associates an already defined Epetra_Vector as the matrix (column ordered) of the linear system.
Referenced by Factor(), and FactorSolve().
int Pliris::SetMatrix | ( | Epetra_SerialDenseVector * | A | ) |
Pliris Matrix Set.
Associates an already defined Epetra_SerialDenseVector as the matrix (column ordered) of the linear system.
int Pliris::SetRHS | ( | Epetra_MultiVector * | B | ) |
Pliris RHS Set.
Associates an already defined Epetra_MultiVector (or Epetra_Vector) as the right-hand-side of the linear system.
int Pliris::Solve | ( | int * | permute, |
int * | num_rhs | ||
) |
Pliris Solve.
Solves the previously factored dense matrix
permute(In) | – permutation matrix |
num_rhs(In) | – factor and solve time in seconds \ Note that the matrix has been previously factored by Factor \ The RHS has been set by SetRHS(Epetra_MultiVector * B) \ On output the result is in the RHS |