EpetraExt  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
List of all members
Epetra_PETScAIJMatrix Class Reference

Epetra_PETScAIJMatrix: A class for constructing and using real-valued sparse compressed row matrices. More...

#include <EpetraExt_PETScAIJMatrix.h>

Inheritance diagram for Epetra_PETScAIJMatrix:
Inheritance graph
[legend]

Constructors/Destructor

 Epetra_PETScAIJMatrix (Mat Amat)
 Epetra_PETScAIJMatrix constructor. More...
 
virtual ~Epetra_PETScAIJMatrix ()
 Epetra_PETScAIJMatrix Destructor. More...
 

Extraction methods

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...
 
int ExtractDiagonalCopy (Epetra_Vector &Diagonal) const
 Returns a copy of the main diagonal in a user-provided vector. More...
 

Computational methods

int Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_PETScAIJMatrix multiplied by a Epetra_MultiVector X in Y. More...
 
int Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_PETScAIJMatrix multiplied by a Epetra_MultiVector X in Y. More...
 
int InvRowSums (Epetra_Vector &x) const
 Computes the sum of absolute values of the rows of the Epetra_PETScAIJMatrix, results returned in x. More...
 
int LeftScale (const Epetra_Vector &x)
 Scales the Epetra_PETScAIJMatrix on the left with a Epetra_Vector x. More...
 
int InvColSums (Epetra_Vector &x) const
 Computes the sum of absolute values of the columns of the Epetra_PETScAIJMatrix, results returned in x. More...
 
int RightScale (const Epetra_Vector &x)
 Scales the Epetra_PETScAIJMatrix on the right with a Epetra_Vector x. More...
 

Matrix Properties Query Methods

bool Filled () const
 If FillComplete() has been called, this query returns true, otherwise it returns false. More...
 
bool LowerTriangular () const
 If matrix is lower triangular, this query returns true, otherwise it returns false. More...
 
bool UpperTriangular () const
 If matrix is upper triangular, this query returns true, otherwise it returns false. More...
 

Attribute access functions

Mat Amat () const
 Returns a pointer to the PETSc matrix used to create this object. More...
 
double NormInf () const
 Returns the infinity norm of the global matrix. More...
 
double NormOne () const
 Returns the one norm of the global matrix. More...
 
int NumGlobalNonzeros () const
 Returns the number of nonzero entries in the global matrix. More...
 
int NumGlobalRows () const
 Returns the number of global matrix rows. More...
 
int NumGlobalCols () const
 Returns the number of global matrix columns. More...
 
int NumGlobalDiagonals () const
 Returns the number of global nonzero diagonal entries. More...
 
long long NumGlobalNonzeros64 () const
 Returns the number of nonzero entries in the global matrix. More...
 
long long NumGlobalRows64 () const
 Returns the number of global matrix rows. More...
 
long long NumGlobalCols64 () const
 Returns the number of global matrix columns. More...
 
long long NumGlobalDiagonals64 () const
 Returns the number of global nonzero diagonal entries. More...
 
int NumMyNonzeros () const
 Returns the number of nonzero entries in the calling processor's portion of the matrix. More...
 
int NumMyRows () const
 Returns the number of matrix rows owned by the calling processor. More...
 
int NumMyCols () const
 Returns the number of matrix columns owned by the calling processor. More...
 
int NumMyDiagonals () const
 Returns the number of local nonzero diagonal entries. More...
 
const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this operator. More...
 
const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this operator (same as domain). More...
 
const Epetra_BlockMapMap () const
 Implement the Epetra_SrcDistObjec::Map() function. More...
 
const Epetra_MapRowMatrixRowMap () const
 Returns the Row Map object needed for implementing Epetra_RowMatrix. More...
 
const Epetra_MapRowMatrixColMap () const
 Returns the Column Map object needed for implementing Epetra_RowMatrix. More...
 
virtual const Epetra_ImportRowMatrixImporter () const
 Returns the Epetra_Import object that contains the import operations for distributed operations. More...
 
const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this matrix. More...
 

I/O Methods

virtual void Print (std::ostream &os) const
 Print method. More...
 

Additional methods required to support the Epetra_Operator interface

const char * Label () const
 Returns a character string describing the operator. More...
 
int SetUseTranspose (bool UseTranspose)
 If set true, transpose of this operator will be applied. More...
 
int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y. More...
 
int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y. More...
 
virtual bool HasNormInf () const
 Returns true because this class can compute an Inf-norm. More...
 
virtual bool UseTranspose () const
 Returns the current UseTranspose setting. More...
 

Additional methods required to implement RowMatrix interface

int NumMyRowEntries (int MyRow, int &NumEntries) const
 Return the current number of values stored for the specified local row. More...
 
int MaxNumEntries () const
 Returns the maximum of NumMyRowEntries() over all rows. More...
 

Additional Inherited Members

Detailed Description

Epetra_PETScAIJMatrix: A class for constructing and using real-valued sparse compressed row matrices.

The Epetra_PETScAIJMatrix is a wrapper class for PETSc sequential or parallel AIJ matrices. It is derived from the Epetra_RowMatrix class, and so provides PETSc users access to Trilinos preconditioners. This class is lightweight, i.e., there are no deep copies of matrix data. Whenever possible, class methods utilize callbacks to native PETSc functions. Currently, only sequential and parallel point AIJ PETSc matrix types are supported.

Definition at line 81 of file EpetraExt_PETScAIJMatrix.h.

Constructor & Destructor Documentation

Epetra_PETScAIJMatrix::Epetra_PETScAIJMatrix ( Mat  Amat)

Epetra_PETScAIJMatrix constructor.

Creates a Epetra_PETScAIJMatrix object by encapsulating an existing PETSc matrix.

Parameters
InAmat - A completely constructed PETSc SEQAIJ or MPIAIJ matrix.
virtual Epetra_PETScAIJMatrix::~Epetra_PETScAIJMatrix ( )
virtual

Epetra_PETScAIJMatrix Destructor.

Member Function Documentation

int Epetra_PETScAIJMatrix::ExtractMyRowCopy ( int  MyRow,
int  Length,
int &  NumEntries,
double *  Values,
int *  Indices 
) const
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.

Implements Epetra_RowMatrix.

int Epetra_PETScAIJMatrix::ExtractDiagonalCopy ( Epetra_Vector Diagonal) const
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.

Implements Epetra_RowMatrix.

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

Returns the result of a Epetra_PETScAIJMatrix 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.

Implements Epetra_RowMatrix.

int Epetra_PETScAIJMatrix::Solve ( bool  Upper,
bool  Trans,
bool  UnitDiagonal,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
virtual

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

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.

Implements Epetra_RowMatrix.

int Epetra_PETScAIJMatrix::InvRowSums ( Epetra_Vector x) const
virtual

Computes the sum of absolute values of the rows of the Epetra_PETScAIJMatrix, 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.

Implements Epetra_RowMatrix.

int Epetra_PETScAIJMatrix::LeftScale ( const Epetra_Vector x)
virtual

Scales the Epetra_PETScAIJMatrix 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.

Implements Epetra_RowMatrix.

int Epetra_PETScAIJMatrix::InvColSums ( Epetra_Vector x) const
virtual

Computes the sum of absolute values of the columns of the Epetra_PETScAIJMatrix, 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.

Implements Epetra_RowMatrix.

int Epetra_PETScAIJMatrix::RightScale ( const Epetra_Vector x)
virtual

Scales the Epetra_PETScAIJMatrix 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.

Implements Epetra_RowMatrix.

bool Epetra_PETScAIJMatrix::Filled ( ) const
inlinevirtual

If FillComplete() has been called, this query returns true, otherwise it returns false.

Implements Epetra_RowMatrix.

Definition at line 214 of file EpetraExt_PETScAIJMatrix.h.

bool Epetra_PETScAIJMatrix::LowerTriangular ( ) const
inlinevirtual

If matrix is lower triangular, this query returns true, otherwise it returns false.

Implements Epetra_RowMatrix.

Definition at line 217 of file EpetraExt_PETScAIJMatrix.h.

bool Epetra_PETScAIJMatrix::UpperTriangular ( ) const
inlinevirtual

If matrix is upper triangular, this query returns true, otherwise it returns false.

Implements Epetra_RowMatrix.

Definition at line 220 of file EpetraExt_PETScAIJMatrix.h.

Mat Epetra_PETScAIJMatrix::Amat ( ) const
inline

Returns a pointer to the PETSc matrix used to create this object.

Definition at line 228 of file EpetraExt_PETScAIJMatrix.h.

double Epetra_PETScAIJMatrix::NormInf ( ) const
virtual

Returns the infinity norm of the global matrix.

Implements Epetra_RowMatrix.

double Epetra_PETScAIJMatrix::NormOne ( ) const
virtual

Returns the one norm of the global matrix.

Implements Epetra_RowMatrix.

int Epetra_PETScAIJMatrix::NumGlobalNonzeros ( ) const
inlinevirtual

Returns the number of nonzero entries in the global matrix.

Implements Epetra_RowMatrix.

Definition at line 244 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::NumGlobalRows ( ) const
inlinevirtual

Returns the number of global matrix rows.

Implements Epetra_RowMatrix.

Definition at line 247 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::NumGlobalCols ( ) const
inlinevirtual

Returns the number of global matrix columns.

Implements Epetra_RowMatrix.

Definition at line 250 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::NumGlobalDiagonals ( ) const
inlinevirtual

Returns the number of global nonzero diagonal entries.

Implements Epetra_RowMatrix.

Definition at line 253 of file EpetraExt_PETScAIJMatrix.h.

long long Epetra_PETScAIJMatrix::NumGlobalNonzeros64 ( ) const
inline

Returns the number of nonzero entries in the global matrix.

Definition at line 256 of file EpetraExt_PETScAIJMatrix.h.

long long Epetra_PETScAIJMatrix::NumGlobalRows64 ( ) const
inline

Returns the number of global matrix rows.

Definition at line 259 of file EpetraExt_PETScAIJMatrix.h.

long long Epetra_PETScAIJMatrix::NumGlobalCols64 ( ) const
inline

Returns the number of global matrix columns.

Definition at line 262 of file EpetraExt_PETScAIJMatrix.h.

long long Epetra_PETScAIJMatrix::NumGlobalDiagonals64 ( ) const
inline

Returns the number of global nonzero diagonal entries.

Definition at line 265 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::NumMyNonzeros ( ) const
inlinevirtual

Returns the number of nonzero entries in the calling processor's portion of the matrix.

Implements Epetra_RowMatrix.

Definition at line 268 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::NumMyRows ( ) const
inlinevirtual

Returns the number of matrix rows owned by the calling processor.

Implements Epetra_RowMatrix.

Definition at line 271 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::NumMyCols ( ) const
inlinevirtual

Returns the number of matrix columns owned by the calling processor.

Implements Epetra_RowMatrix.

Definition at line 274 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::NumMyDiagonals ( ) const
inlinevirtual

Returns the number of local nonzero diagonal entries.

Implements Epetra_RowMatrix.

Definition at line 277 of file EpetraExt_PETScAIJMatrix.h.

const Epetra_Map& Epetra_PETScAIJMatrix::OperatorDomainMap ( ) const
inlinevirtual

Returns the Epetra_Map object associated with the domain of this operator.

Implements Epetra_Operator.

Definition at line 280 of file EpetraExt_PETScAIJMatrix.h.

const Epetra_Map& Epetra_PETScAIJMatrix::OperatorRangeMap ( ) const
inlinevirtual

Returns the Epetra_Map object associated with the range of this operator (same as domain).

Implements Epetra_Operator.

Definition at line 283 of file EpetraExt_PETScAIJMatrix.h.

const Epetra_BlockMap& Epetra_PETScAIJMatrix::Map ( ) const
inlinevirtual

Implement the Epetra_SrcDistObjec::Map() function.

Implements Epetra_RowMatrix.

Definition at line 286 of file EpetraExt_PETScAIJMatrix.h.

const Epetra_Map& Epetra_PETScAIJMatrix::RowMatrixRowMap ( ) const
inlinevirtual

Returns the Row Map object needed for implementing Epetra_RowMatrix.

Implements Epetra_RowMatrix.

Definition at line 289 of file EpetraExt_PETScAIJMatrix.h.

const Epetra_Map& Epetra_PETScAIJMatrix::RowMatrixColMap ( ) const
inlinevirtual

Returns the Column Map object needed for implementing Epetra_RowMatrix.

Implements Epetra_RowMatrix.

Definition at line 292 of file EpetraExt_PETScAIJMatrix.h.

virtual const Epetra_Import* Epetra_PETScAIJMatrix::RowMatrixImporter ( ) const
inlinevirtual

Returns the Epetra_Import object that contains the import operations for distributed operations.

Implements Epetra_RowMatrix.

Definition at line 295 of file EpetraExt_PETScAIJMatrix.h.

const Epetra_Comm& Epetra_PETScAIJMatrix::Comm ( ) const
inlinevirtual

Returns a pointer to the Epetra_Comm communicator associated with this matrix.

Implements Epetra_Operator.

Definition at line 298 of file EpetraExt_PETScAIJMatrix.h.

virtual void Epetra_PETScAIJMatrix::Print ( std::ostream &  os) const
virtual

Print method.

const char* Epetra_PETScAIJMatrix::Label ( ) const
inlinevirtual

Returns a character string describing the operator.

Implements Epetra_Operator.

Definition at line 313 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::SetUseTranspose ( bool  UseTranspose)
inlinevirtual

If set true, transpose of this operator will be applied.

  This flag allows the transpose of the given operator to be used implicitly.  Setting this flag
  affects only the Apply() and ApplyInverse() methods.  If the implementation of this interface 

does not support transpose use, this method should return a value of -1.

Parameters
InUseTranspose -If true, multiply by the transpose of operator, otherwise just use operator.
Returns
Always returns 0.

Implements Epetra_Operator.

Definition at line 325 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::Apply ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.

Parameters
InX - A Epetra_MultiVector of dimension NumVectors to multiply with matrix.
OutY -A Epetra_MultiVector of dimension NumVectors containing result.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_Operator.

Definition at line 337 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::ApplyInverse ( const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

Returns the result of a Epetra_Operator inverse applied to an Epetra_MultiVector X in Y.

  In this implementation, we use several existing attributes to determine how virtual
  method ApplyInverse() should call the concrete method Solve().  We pass in the UpperTriangular(), 

the Epetra_PETScAIJMatrix::UseTranspose(), and NoDiagonal() methods. The most notable warning is that if a matrix has no diagonal values we assume that there is an implicit unit diagonal that should be accounted for when doing a triangular solve.

Parameters
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.

Implements Epetra_Operator.

Definition at line 354 of file EpetraExt_PETScAIJMatrix.h.

virtual bool Epetra_PETScAIJMatrix::HasNormInf ( ) const
inlinevirtual

Returns true because this class can compute an Inf-norm.

Implements Epetra_Operator.

Definition at line 359 of file EpetraExt_PETScAIJMatrix.h.

virtual bool Epetra_PETScAIJMatrix::UseTranspose ( ) const
inlinevirtual

Returns the current UseTranspose setting.

Implements Epetra_Operator.

Definition at line 362 of file EpetraExt_PETScAIJMatrix.h.

int Epetra_PETScAIJMatrix::NumMyRowEntries ( int  MyRow,
int &  NumEntries 
) const
virtual

Return the current number of values stored for the specified local row.

Similar to NumMyEntries() except NumEntries is returned as an argument
and error checking is done on the input value MyRow.
Parameters
InMyRow - Local row.
OutNumEntries - Number of nonzero values.
Returns
Integer error code, set to 0 if successful.

Implements Epetra_RowMatrix.

int Epetra_PETScAIJMatrix::MaxNumEntries ( ) const
virtual

Returns the maximum of NumMyRowEntries() over all rows.

Implements Epetra_RowMatrix.


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