Epetra  Development
 All Classes Namespaces Files Functions Variables Enumerations Enumerator Pages
List of all members
Epetra_RowMatrix Class Referenceabstract

Epetra_RowMatrix: A pure virtual class for using real-valued double-precision row matrices. More...

#include <Epetra_RowMatrix.h>

Inheritance diagram for Epetra_RowMatrix:
Inheritance graph
[legend]
Collaboration diagram for Epetra_RowMatrix:
Collaboration graph
[legend]

Public Member Functions

Destructor
virtual ~Epetra_RowMatrix ()
 Destructor.
 
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.
 
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.
 
virtual double NormInf () const =0
 Returns the infinity norm of the global matrix.
 
virtual double NormOne () const =0
 Returns the one norm of the global matrix.
 
virtual int NumGlobalNonzeros () const =0
 Returns the number of nonzero entries in the global matrix.
 
virtual long long NumGlobalNonzeros64 () const =0
 
virtual int NumGlobalRows () const =0
 Returns the number of global matrix rows.
 
virtual long long NumGlobalRows64 () const =0
 
virtual int NumGlobalCols () const =0
 Returns the number of global matrix columns.
 
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.
 
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.
 
virtual int NumMyRows () const =0
 Returns the number of matrix rows owned by the calling processor.
 
virtual int NumMyCols () const =0
 Returns the number of matrix columns owned by the calling processor.
 
virtual int NumMyDiagonals () const =0
 Returns the number of local nonzero diagonal entries, based on global row/column index comparisons.
 
virtual bool LowerTriangular () const =0
 If matrix is lower triangular in local index space, this query returns true, otherwise it returns false.
 
virtual bool UpperTriangular () const =0
 If matrix is upper triangular in local index space, this query returns true, otherwise it returns false.
 
virtual const Epetra_MapRowMatrixRowMap () const =0
 Returns the Epetra_Map object associated with the rows of this matrix.
 
virtual const Epetra_MapRowMatrixColMap () const =0
 Returns the Epetra_Map object associated with the columns of this matrix.
 
virtual const Epetra_ImportRowMatrixImporter () const =0
 Returns the Epetra_Import object that contains the import operations for distributed operations.
 
- Public Member Functions inherited from Epetra_Operator
virtual ~Epetra_Operator ()
 Destructor.
 
virtual int SetUseTranspose (bool UseTranspose)=0
 If set true, transpose of this operator will be applied. More...
 
virtual int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
 Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y. More...
 
virtual int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const =0
 Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y. More...
 
virtual const char * Label () const =0
 Returns a character string describing the operator.
 
virtual bool UseTranspose () const =0
 Returns the current UseTranspose setting.
 
virtual bool HasNormInf () const =0
 Returns true if the this object can provide an approximate Inf-norm, false otherwise.
 
virtual const Epetra_CommComm () const =0
 Returns a pointer to the Epetra_Comm communicator associated with this operator.
 
virtual const Epetra_MapOperatorDomainMap () const =0
 Returns the Epetra_Map object associated with the domain of this operator.
 
virtual const Epetra_MapOperatorRangeMap () const =0
 Returns the Epetra_Map object associated with the range of this operator.
 
- Public Member Functions inherited from Epetra_SrcDistObject
virtual ~Epetra_SrcDistObject ()
 Epetra_SrcDistObject destructor.
 
virtual const Epetra_BlockMapMap () const =0
 Returns a reference to the Epetra_BlockMap for this object.
 

Detailed Description

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.

Member Function Documentation

virtual int Epetra_RowMatrix::ExtractDiagonalCopy ( Epetra_Vector Diagonal) const
pure virtual

Returns a copy of the main diagonal in a user-provided vector.

Parameters
OutDiagonal - Extracted main diagonal.
Returns
Integer error code, set to 0 if successful.

Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_BasicRowMatrix, Epetra_OskiMatrix, and Epetra_MsrMatrix.

virtual int Epetra_RowMatrix::ExtractMyRowCopy ( int  MyRow,
int  Length,
int &  NumEntries,
double *  Values,
int *  Indices 
) const
pure virtual

Returns a copy of the specified local row in user-provided arrays.

Parameters
InMyRow - Local row to extract.
InLength - Length of Values and Indices.
OutNumEntries - Number of nonzero entries extracted.
OutValues - Extracted values for this row.
OutIndices - Extracted local column indices for the corresponding values.
Returns
Integer error code, set to 0 if successful.

Implemented in Epetra_VbrMatrix, Epetra_CrsMatrix, Epetra_BasicRowMatrix, Epetra_VbrRowMatrix, Epetra_JadMatrix, and Epetra_MsrMatrix.

virtual int Epetra_RowMatrix::InvColSums ( Epetra_Vector x) const
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.

Parameters
Outx -A Epetra_Vector containing the column sums of the this matrix.
Warning
It is assumed that the distribution of x is the same as the rows of this.
Returns
Integer error code, set to 0 if successful.

Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_BasicRowMatrix, and Epetra_MsrMatrix.

virtual int Epetra_RowMatrix::InvRowSums ( Epetra_Vector x) const
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.

Parameters
Outx -A Epetra_Vector containing the row sums of the this matrix.
Warning
It is assumed that the distribution of x is the same as the rows of this.
Returns
Integer error code, set to 0 if successful.

Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_BasicRowMatrix, and Epetra_MsrMatrix.

virtual int Epetra_RowMatrix::LeftScale ( const Epetra_Vector x)
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.
Parameters
Inx -A Epetra_Vector to solve for.
Returns
Integer error code, set to 0 if successful.

Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_BasicRowMatrix, Epetra_VbrRowMatrix, and Epetra_MsrMatrix.

virtual int Epetra_RowMatrix::Multiply ( bool  TransA,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
pure virtual

Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y.

Parameters
InTransA -If true, multiply by the transpose of matrix, otherwise just use matrix.
InX - A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
OutY -A Epetra_MultiVector of dimension NumVectorscontaining result.
Returns
Integer error code, set to 0 if successful.

Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_OskiMatrix, Epetra_BasicRowMatrix, Epetra_VbrRowMatrix, Epetra_JadMatrix, and Epetra_MsrMatrix.

Referenced by Epetra_BasicRowMatrix::Apply().

virtual int Epetra_RowMatrix::NumMyRowEntries ( int  MyRow,
int &  NumEntries 
) const
pure virtual

Returns the number of nonzero entries in MyRow.

Parameters
InMyRow - Local row.
OutNumEntries - Number of nonzero values present.
Returns
Integer error code, set to 0 if successful.

Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_MsrMatrix, Epetra_BasicRowMatrix, Epetra_VbrRowMatrix, and Epetra_JadMatrix.

virtual int Epetra_RowMatrix::RightScale ( const Epetra_Vector x)
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.
Parameters
Inx -The Epetra_Vector used for scaling this.
Returns
Integer error code, set to 0 if successful.

Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_BasicRowMatrix, Epetra_MsrMatrix, and Epetra_VbrRowMatrix.

virtual int Epetra_RowMatrix::Solve ( bool  Upper,
bool  Trans,
bool  UnitDiagonal,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
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.
Parameters
InUpper -If true, solve Ux = y, otherwise solve Lx = y.
InTrans -If true, solve transpose problem.
InUnitDiagonal -If true, assume diagonal is unit (whether it's stored or not).
InX - A Epetra_MultiVector of dimension NumVectors to solve for.
OutY -A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Implemented in Epetra_CrsMatrix, Epetra_VbrMatrix, Epetra_OskiMatrix, Epetra_BasicRowMatrix, Epetra_VbrRowMatrix, Epetra_JadMatrix, and Epetra_MsrMatrix.

Referenced by Epetra_CrsMatrix::ApplyInverse().


The documentation for this class was generated from the following file: