Ifpack_LocalFilter a class for light-weight extraction of the submatrix corresponding to local rows and columns. More...
#include <Ifpack_LocalFilter.h>
Public Member Functions | |
int | SetOwnership (bool) |
Sets ownership. | |
int | SetUseTranspose (bool UseTranspose_in) |
Sets use transpose (not implemented). | |
bool | UseTranspose () const |
Returns the current UseTranspose setting. | |
bool | HasNormInf () const |
Returns true if the this object can provide an approximate Inf-norm, false otherwise. | |
const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this operator. | |
const Epetra_Map & | OperatorDomainMap () const |
Returns the Epetra_Map object associated with the domain of this operator. | |
const Epetra_Map & | OperatorRangeMap () const |
Returns the Epetra_Map object associated with the range of this operator. | |
const Epetra_BlockMap & | Map () const |
const char * | Label () const |
Ifpack_LocalFilter (const Teuchos::RefCountPtr< const Epetra_RowMatrix > &Matrix) | |
Constructor. | |
virtual | ~Ifpack_LocalFilter () |
Destructor. | |
virtual int | NumMyRowEntries (int MyRow, int &NumEntries) const |
Returns the number of nonzero entries in MyRow. More... | |
virtual int | MaxNumEntries () const |
Returns the maximum of NumMyRowEntries() over all rows. | |
virtual 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... | |
virtual int | ExtractDiagonalCopy (Epetra_Vector &Diagonal) const |
Returns a copy of the main diagonal in a user-provided vector. More... | |
virtual int | Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y. More... | |
virtual int | Solve (bool, bool, bool, const Epetra_MultiVector &, Epetra_MultiVector &) const |
Returns result of a local-only solve using a triangular Epetra_RowMatrix with Epetra_MultiVectors X and Y (NOT IMPLEMENTED). | |
virtual int | Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
virtual int | ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
virtual int | InvRowSums (Epetra_Vector &) const |
Computes the sum of absolute values of the rows of the Epetra_RowMatrix, results returned in x (NOT IMPLEMENTED). | |
virtual int | LeftScale (const Epetra_Vector &) |
Scales the Epetra_RowMatrix on the left with a Epetra_Vector x (NOT IMPLEMENTED). | |
virtual int | InvColSums (Epetra_Vector &) const |
Computes the sum of absolute values of the columns of the Epetra_RowMatrix, results returned in x (NOT IMPLEMENTED). | |
virtual int | RightScale (const Epetra_Vector &) |
Scales the Epetra_RowMatrix on the right with a Epetra_Vector x (NOT IMPLEMENTED). | |
virtual bool | Filled () const |
If FillComplete() has been called, this query returns true, otherwise it returns false. | |
virtual double | NormInf () const |
Returns the infinity norm of the global matrix. | |
virtual double | NormOne () const |
Returns the one norm of the global matrix. | |
virtual int | NumGlobalNonzeros () const |
Returns the number of nonzero entries in the global matrix. | |
virtual int | NumGlobalRows () const |
Returns the number of global matrix rows. | |
virtual int | NumGlobalCols () const |
Returns the number of global matrix columns. | |
virtual int | NumGlobalDiagonals () const |
Returns the number of global nonzero diagonal entries, based on global row/column index comparisons. | |
virtual long long | NumGlobalNonzeros64 () const |
Returns the number of nonzero entries in the global matrix. | |
virtual long long | NumGlobalRows64 () const |
Returns the number of global matrix rows. | |
virtual long long | NumGlobalCols64 () const |
Returns the number of global matrix columns. | |
virtual long long | NumGlobalDiagonals64 () const |
Returns the number of global nonzero diagonal entries, based on global row/column index comparisons. | |
virtual int | NumMyNonzeros () const |
Returns the number of nonzero entries in the calling processor's portion of the matrix. | |
virtual int | NumMyRows () const |
Returns the number of matrix rows owned by the calling processor. | |
virtual int | NumMyCols () const |
Returns the number of matrix columns owned by the calling processor. | |
virtual int | NumMyDiagonals () const |
Returns the number of local nonzero diagonal entries, based on global row/column index comparisons. | |
virtual bool | LowerTriangular () const |
If matrix is lower triangular in local index space, this query returns true, otherwise it returns false. | |
virtual bool | UpperTriangular () const |
If matrix is upper triangular in local index space, this query returns true, otherwise it returns false. | |
virtual const Epetra_Map & | RowMatrixRowMap () const |
Returns the Epetra_Map object associated with the rows of this matrix. | |
virtual const Epetra_Map & | RowMatrixColMap () const |
Returns the Epetra_Map object associated with the columns of this matrix. | |
virtual const Epetra_Import * | RowMatrixImporter () const |
Returns the Epetra_Import object that contains the import operations for distributed operations. | |
Ifpack_LocalFilter a class for light-weight extraction of the submatrix corresponding to local rows and columns.
Class Ifpack_LocalFilter enables a light-weight contruction of an Epetra_RowMatrix-derived object, containing only the elements of the original, distributed matrix with local row and column ID. The local submatrix is based on a communicator containing the local process only. Each process will have its local object, corresponding to the local submatrix. Submatrices may or may not overlap.
The following instructions can be used to create "localized" matrices:
Once created, LocalA
defined, on each process, the submatrix corresponding to local rows and columns only. The creation and use of LocalA
is "cheap", as the elements of the local matrix are obtained through calls to ExtractMyRowCopy on the original, distributed matrix, say A. This means that A
must remain in scope every time LocalA
is accessed.
A very convenient use of this class is to use Ifpack solvers to compute the LU factorizations of local blocks. If applied to a localized matrix, several Ifpack objects can operator in the same phase in a safe way, without non-required data exchange.
Definition at line 98 of file Ifpack_LocalFilter.h.
|
virtual |
Returns a copy of the main diagonal in a user-provided vector.
Diagonal | - (Out) Extracted main diagonal. |
Implements Epetra_RowMatrix.
Definition at line 167 of file Ifpack_LocalFilter.cpp.
|
virtual |
Returns a copy of the specified local row in user-provided arrays.
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 global column indices for the corresponding values. |
Implements Epetra_RowMatrix.
Definition at line 131 of file Ifpack_LocalFilter.cpp.
Referenced by Ifpack_LocalFilter().
|
inlinevirtual |
Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y.
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. |
Implements Epetra_RowMatrix.
Definition at line 175 of file Ifpack_LocalFilter.h.
|
inlinevirtual |
Returns the number of nonzero entries in MyRow.
MyRow | - (In) Local row. |
NumEntries | - (Out) Number of nonzero values present. |
Implements Epetra_RowMatrix.
Definition at line 123 of file Ifpack_LocalFilter.h.