Ifpack Package Browser (Single Doxygen Collection)
Development
|
Ifpack_Container: a pure virtual class for creating and solving local linear problems. More...
#include <Ifpack_Container.h>
Public Member Functions | |
virtual | ~Ifpack_Container () |
Destructor. More... | |
virtual int | NumRows () const =0 |
Returns the number of rows of the matrix and LHS/RHS. More... | |
virtual int | NumVectors () const =0 |
Returns the number of vectors in LHS/RHS. More... | |
virtual int | SetNumVectors (const int i)=0 |
Sets the number of vectors for LHS/RHS. More... | |
virtual double & | LHS (const int i, const int Vector=0)=0 |
Returns the i-th component of the vector Vector of LHS. More... | |
virtual double & | RHS (const int i, const int Vector=0)=0 |
Returns the i-th component of the vector Vector of RHS. More... | |
virtual int & | ID (const int i)=0 |
Returns the ID associated to local row i. More... | |
virtual int | SetMatrixElement (const int row, const int col, const double value)=0 |
Set the matrix element (row,col) to value . More... | |
virtual int | Initialize ()=0 |
Initializes the container, by performing all operations that only require matrix structure. More... | |
virtual int | Compute (const Epetra_RowMatrix &A)=0 |
Finalizes the linear system matrix and prepares for the application of the inverse. More... | |
virtual int | SetParameters (Teuchos::ParameterList &List)=0 |
Sets all necessary parameters. More... | |
virtual bool | IsInitialized () const =0 |
Returns true is the container has been successfully initialized. More... | |
virtual bool | IsComputed () const =0 |
Returns true is the container has been successfully computed. More... | |
virtual int | Apply ()=0 |
Apply the matrix to RHS, results are stored in LHS. More... | |
virtual int | ApplyInverse ()=0 |
Apply the inverse of the matrix to RHS, results are stored in LHS. More... | |
virtual const char * | Label () const =0 |
Returns the label of this container. More... | |
virtual double | InitializeFlops () const =0 |
Returns the flops in Initialize(). More... | |
virtual double | ComputeFlops () const =0 |
Returns the flops in Compute(). More... | |
virtual double | ApplyFlops () const =0 |
Returns the flops in Apply(). More... | |
virtual double | ApplyInverseFlops () const =0 |
Returns the flops in ApplyInverse(). More... | |
virtual std::ostream & | Print (std::ostream &os) const =0 |
Prints out basic information about the container. More... | |
Ifpack_Container: a pure virtual class for creating and solving local linear problems.
Class Ifpack_Container provides the abstract interfaces for containers. A "container" is an object that hosts all it is necessary to create, populate, and solve local linear problems. The local linear problem matrix, B, is a submatrix of the local components of a distributed matrix, A. The idea of container is to specify the rows of A that are contained in B, then extract B from A, and compute all it is necessary to solve a linear system in B. Then, set starting solution (if necessary) and right-hand side for B, and solve the linear system in B.
A container should be used in the following manner:
The number of vectors can be set using SetNumVectors(), and it is defaulted to 1.
Containers are currently used by class Ifpack_BlockRelaxation.
Ifpack_Container is a pure virtual class. Two concrete implementations are provided in classes Ifpack_SparseContainer (that stores matrices in sparse the format Epetra_CrsMatrix) and Ifpack_DenseContainer (for relatively small matrices, as matrices are stored as Epetra_SerialDenseMatrix's).
Definition at line 98 of file Ifpack_Container.h.
|
inlinevirtual |
Destructor.
Definition at line 103 of file Ifpack_Container.h.
|
pure virtual |
Returns the number of rows of the matrix and LHS/RHS.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns the number of vectors in LHS/RHS.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Sets the number of vectors for LHS/RHS.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns the i-th component of the vector Vector of LHS.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns the i-th component of the vector Vector of RHS.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns the ID associated to local row i.
The set of (local) rows assigned to this container is defined by calling ID(i) = j, where i (from 0 to NumRows()) indicates the container-row, and j indicates the local row in the calling process.
This is usually used to recorder the local row ID (on calling process) of the i-th row in the container.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Set the matrix element (row,col) to value
.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Initializes the container, by performing all operations that only require matrix structure.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Finalizes the linear system matrix and prepares for the application of the inverse.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Sets all necessary parameters.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns true
is the container has been successfully initialized.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns true
is the container has been successfully computed.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Apply the matrix to RHS, results are stored in LHS.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Apply the inverse of the matrix to RHS, results are stored in LHS.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns the label of this container.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns the flops in Initialize().
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns the flops in Compute().
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns the flops in Apply().
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Returns the flops in ApplyInverse().
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.
|
pure virtual |
Prints out basic information about the container.
Implemented in Ifpack_TriDiContainer, Ifpack_DenseContainer, and Ifpack_SparseContainer< T >.