| Epetra Package Browser (Single Doxygen Collection)
    Development
    | 
Epetra_RowMatrix: A pure virtual class for using real-valued double-precision row matrices. More...
#include <Epetra_RowMatrix.h>

| Destructor | |
| virtual | ~Epetra_RowMatrix () | 
| Destructor.  More... | |
| Matrix data extraction routines | |
| virtual int | NumMyRowEntries (int MyRow, int &NumEntries) const =0 | 
| Returns the number of nonzero entries in MyRow.  More... | |
| virtual int | MaxNumEntries () const =0 | 
| Returns the maximum of NumMyRowEntries() over all rows.  More... | |
| virtual int | ExtractMyRowCopy (int MyRow, int Length, int &NumEntries, double *Values, int *Indices) const =0 | 
| Returns a copy of the specified local row in user-provided arrays.  More... | |
| virtual int | ExtractDiagonalCopy (Epetra_Vector &Diagonal) const =0 | 
| Returns a copy of the main diagonal in a user-provided vector.  More... | |
| Mathematical functions | |
| virtual int | Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0 | 
| Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y.  More... | |
| virtual int | Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0 | 
| Returns result of a local-only solve using a triangular Epetra_RowMatrix with Epetra_MultiVectors X and Y.  More... | |
| virtual int | InvRowSums (Epetra_Vector &x) const =0 | 
| Computes the sum of absolute values of the rows of the Epetra_RowMatrix, results returned in x.  More... | |
| virtual int | LeftScale (const Epetra_Vector &x)=0 | 
| Scales the Epetra_RowMatrix on the left with a Epetra_Vector x.  More... | |
| virtual int | InvColSums (Epetra_Vector &x) const =0 | 
| Computes the sum of absolute values of the columns of the Epetra_RowMatrix, results returned in x.  More... | |
| virtual int | RightScale (const Epetra_Vector &x)=0 | 
| Scales the Epetra_RowMatrix on the right with a Epetra_Vector x.  More... | |
| Attribute access functions | |
| virtual bool | Filled () const =0 | 
| If FillComplete() has been called, this query returns true, otherwise it returns false.  More... | |
| virtual double | NormInf () const =0 | 
| Returns the infinity norm of the global matrix.  More... | |
| virtual double | NormOne () const =0 | 
| Returns the one norm of the global matrix.  More... | |
| virtual int | NumGlobalNonzeros () const =0 | 
| Returns the number of nonzero entries in the global matrix.  More... | |
| virtual long long | NumGlobalNonzeros64 () const =0 | 
| virtual int | NumGlobalRows () const =0 | 
| Returns the number of global matrix rows.  More... | |
| virtual long long | NumGlobalRows64 () const =0 | 
| virtual int | NumGlobalCols () const =0 | 
| Returns the number of global matrix columns.  More... | |
| virtual long long | NumGlobalCols64 () const =0 | 
| virtual int | NumGlobalDiagonals () const =0 | 
| Returns the number of global nonzero diagonal entries, based on global row/column index comparisons.  More... | |
| virtual long long | NumGlobalDiagonals64 () const =0 | 
| virtual int | NumMyNonzeros () const =0 | 
| Returns the number of nonzero entries in the calling processor's portion of the matrix.  More... | |
| virtual int | NumMyRows () const =0 | 
| Returns the number of matrix rows owned by the calling processor.  More... | |
| virtual int | NumMyCols () const =0 | 
| Returns the number of matrix columns owned by the calling processor.  More... | |
| virtual int | NumMyDiagonals () const =0 | 
| Returns the number of local nonzero diagonal entries, based on global row/column index comparisons.  More... | |
| virtual bool | LowerTriangular () const =0 | 
| If matrix is lower triangular in local index space, this query returns true, otherwise it returns false.  More... | |
| virtual bool | UpperTriangular () const =0 | 
| If matrix is upper triangular in local index space, this query returns true, otherwise it returns false.  More... | |
| virtual const Epetra_Map & | RowMatrixRowMap () const =0 | 
| Returns the Epetra_Map object associated with the rows of this matrix.  More... | |
| virtual const Epetra_Map & | RowMatrixColMap () const =0 | 
| Returns the Epetra_Map object associated with the columns of this matrix.  More... | |
| virtual const Epetra_Import * | RowMatrixImporter () const =0 | 
| Returns the Epetra_Import object that contains the import operations for distributed operations.  More... | |
Epetra_RowMatrix: A pure virtual class for using real-valued double-precision row matrices.
The Epetra_RowMatrix class is a pure virtual class (specifies interface only) that enable the use of real-valued double-precision sparse matrices where matrix entries are intended for row access. It is currently implemented by both the Epetra_CrsMatrix and Epetra_VbrMatrix classes.
Definition at line 76 of file Epetra_RowMatrix.h.
| 
 | inlinevirtual | 
Destructor.
Definition at line 82 of file Epetra_RowMatrix.h.
| 
 | pure virtual | 
Returns the number of nonzero entries in MyRow.
| In | MyRow - Local row. | 
| Out | NumEntries - Number of nonzero values present. | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_BasicRowMatrix, Epetra_VbrRowMatrix, and Epetra_JadMatrix.
| 
 | pure virtual | 
Returns the maximum of NumMyRowEntries() over all rows.
Implemented in Epetra_VbrMatrix, Epetra_CrsMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns a copy of the specified local row in user-provided arrays.
| In | MyRow - Local row to extract. | 
| In | Length - Length of Values and Indices. | 
| Out | NumEntries - Number of nonzero entries extracted. | 
| Out | Values - Extracted values for this row. | 
| Out | Indices - Extracted local column indices for the corresponding values. | 
Implemented in Epetra_VbrMatrix, Epetra_CrsMatrix, Epetra_BasicRowMatrix, Epetra_VbrRowMatrix, and Epetra_JadMatrix.
| 
 | pure virtual | 
Returns a copy of the main diagonal in a user-provided vector.
| Out | Diagonal - Extracted main diagonal. | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_BasicRowMatrix, and Epetra_OskiMatrix.
| 
 | pure virtual | 
Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y.
| In | TransA -If true, multiply by the transpose of matrix, otherwise just use matrix. | 
| In | X - A Epetra_MultiVector of dimension NumVectors to multiply with matrix. | 
| Out | Y -A Epetra_MultiVector of dimension NumVectorscontaining result. | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_OskiMatrix, Epetra_BasicRowMatrix, Epetra_VbrRowMatrix, and Epetra_JadMatrix.
| 
 | pure virtual | 
Returns result of a local-only solve using a triangular Epetra_RowMatrix with Epetra_MultiVectors X and Y.
This method will perform a triangular solve independently on each processor of the parallel machine. No communication is performed.
| In | Upper -If true, solve Ux = y, otherwise solve Lx = y. | 
| In | Trans -If true, solve transpose problem. | 
| In | UnitDiagonal -If true, assume diagonal is unit (whether it's stored or not). | 
| In | X - A Epetra_MultiVector of dimension NumVectors to solve for. | 
| Out | Y -A Epetra_MultiVector of dimension NumVectors containing result. | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_OskiMatrix, Epetra_BasicRowMatrix, Epetra_VbrRowMatrix, and Epetra_JadMatrix.
| 
 | pure virtual | 
Computes the sum of absolute values of the rows of the Epetra_RowMatrix, results returned in x.
The vector x will return such that x[i] will contain the inverse of sum of the absolute values of the \e this matrix will be scaled such that A(i,j) = x(i)*A(i,j) where i denotes the global row number of A and j denotes the global column number of A. Using the resulting vector from this function as input to LeftScale()
will make the infinity norm of the resulting matrix exactly 1.
| Out | x -A Epetra_Vector containing the row sums of the this matrix. | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Scales the Epetra_RowMatrix on the left with a Epetra_Vector x.
The \e 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.
| In | x -A Epetra_Vector to solve for. | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_BasicRowMatrix, and Epetra_VbrRowMatrix.
| 
 | pure virtual | 
Computes the sum of absolute values of the columns of the Epetra_RowMatrix, results returned in x.
The vector x will return such that x[j] will contain the inverse of sum of the absolute values of the \e this matrix will be sca 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. Using the resulting vector from this function as input to
RighttScale() will make the one norm of the resulting matrix exactly 1.
| Out | x -A Epetra_Vector containing the column sums of the this matrix. | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Scales the Epetra_RowMatrix on the right with a Epetra_Vector x.
The \e 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.
| In | x -The Epetra_Vector used for scaling this. | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_BasicRowMatrix, and Epetra_VbrRowMatrix.
| 
 | pure virtual | 
If FillComplete() has been called, this query returns true, otherwise it returns false.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the infinity norm of the global matrix.
Implements Epetra_Operator.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the one norm of the global matrix.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the number of nonzero entries in the global matrix.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the number of global matrix rows.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the number of global matrix columns.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the number of global nonzero diagonal entries, based on global row/column index comparisons.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the number of nonzero entries in the calling processor's portion of the matrix.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the number of matrix rows owned by the calling processor.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the number of matrix columns owned by the calling processor.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the number of local nonzero diagonal entries, based on global row/column index comparisons.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
If matrix is lower triangular in local index space, this query returns true, otherwise it returns false.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
If matrix is upper triangular in local index space, this query returns true, otherwise it returns false.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the Epetra_Map object associated with the rows of this matrix.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the Epetra_Map object associated with the columns of this matrix.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
| 
 | pure virtual | 
Returns the Epetra_Import object that contains the import operations for distributed operations.
Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, and Epetra_BasicRowMatrix.
 1.8.5
 1.8.5