Epetra Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Private Attributes | List of all members
Epetra_VbrRowMatrix Class Reference

Epetra_VbrRowMatrix: A class for using an existing Epetra_VbrMatrix object as an Epetra_RowMatrix object. More...

#include <Epetra_VbrRowMatrix.h>

Inheritance diagram for Epetra_VbrRowMatrix:
Inheritance graph
[legend]

Private Attributes

Epetra_VbrMatrixmatrix_
 

Constructors/Destructor

 Epetra_VbrRowMatrix (Epetra_VbrMatrix *Matrix)
 Epetra_VbrRowMatrix constuctor. More...
 
virtual ~Epetra_VbrRowMatrix ()
 Epetra_VbrRowMatrix Destructor. More...
 

Post-construction modifications

int UpdateMatrix (Epetra_VbrMatrix *Matrix)
 Update the matrix to which this object points. More...
 

Methods required for implementing Epetra_BasicRowMatrix

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 ExtractMyEntryView (int CurEntry, double *&Value, int &RowIndex, int &ColIndex)
 Returns a reference to the ith entry in the matrix, along with its row and column index. More...
 
int ExtractMyEntryView (int CurEntry, double const *&Value, int &RowIndex, int &ColIndex) const
 Returns a const reference to the ith entry in the matrix, along with its row and column index. More...
 
int NumMyRowEntries (int MyRow, int &NumEntries) const
 Return the current number of values stored for the specified local row. More...
 

Computational methods

int RightScale (const Epetra_Vector &x)
 Scales the Epetra_VbrMatrix on the right with a Epetra_Vector x. More...
 
int LeftScale (const Epetra_Vector &x)
 Scales the Epetra_VbrMatrix on the left with a Epetra_Vector x. More...
 
int Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_VbrRowMatrix 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_VbrRowMatrix solve with a Epetra_MultiVector X in Y (not implemented). More...
 

Additional Inherited Members

- Public Member Functions inherited from Epetra_BasicRowMatrix
 Epetra_BasicRowMatrix (const Epetra_Comm &Comm)
 Epetra_BasicRowMatrix constructor. More...
 
virtual ~Epetra_BasicRowMatrix ()
 Epetra_BasicRowMatrix Destructor. More...
 
void SetMaps (const Epetra_Map &RowMap, const Epetra_Map &ColMap)
 Set maps (Version 1); call this function or the next, but not both. More...
 
void SetMaps (const Epetra_Map &RowMap, const Epetra_Map &ColMap, const Epetra_Map &DomainMap, const Epetra_Map &RangeMap)
 Set maps (Version 2); call this function or the previous, but not both. More...
 
virtual int ExtractDiagonalCopy (Epetra_Vector &Diagonal) const
 Returns a copy of the main diagonal in a user-provided vector. More...
 
virtual int InvRowSums (Epetra_Vector &x) const
 Computes the sum of absolute values of the rows of the Epetra_BasicRowMatrix, results returned in x. More...
 
virtual int InvColSums (Epetra_Vector &x) const
 Computes the sum of absolute values of the columns of the Epetra_BasicRowMatrix, results returned in x. More...
 
virtual bool Filled () const
 If FillComplete() has been called, this query returns true, otherwise it returns false, presently always returns true. More...
 
bool LowerTriangular () const
 If matrix is lower triangular, this query returns true, otherwise it returns false. More...
 
virtual bool UpperTriangular () const
 If matrix is upper triangular, this query returns true, otherwise it returns false. More...
 
virtual double NormInf () const
 Returns the infinity norm of the global matrix. More...
 
virtual double NormOne () const
 Returns the one norm of the global matrix. More...
 
virtual int NumGlobalNonzeros () const
 Returns the number of nonzero entries in the global matrix. More...
 
virtual long long NumGlobalNonzeros64 () const
 
virtual int NumGlobalRows () const
 Returns the number of global matrix rows. More...
 
virtual long long NumGlobalRows64 () const
 
virtual int NumGlobalCols () const
 Returns the number of global matrix columns. More...
 
virtual long long NumGlobalCols64 () const
 
virtual int NumGlobalDiagonals () const
 Returns the number of global nonzero diagonal entries. More...
 
virtual long long NumGlobalDiagonals64 () const
 
virtual int NumMyNonzeros () const
 Returns the number of nonzero entries in the calling processor's portion of the matrix. More...
 
virtual int NumMyRows () const
 Returns the number of matrix rows owned by the calling processor. More...
 
virtual int NumMyCols () const
 Returns the number of matrix columns owned by the calling processor. More...
 
virtual int NumMyDiagonals () const
 Returns the number of local nonzero diagonal entries. More...
 
virtual int MaxNumEntries () const
 Returns the maximum number of nonzero entries across all rows on this processor. More...
 
virtual const Epetra_MapOperatorDomainMap () const
 Returns the Epetra_Map object associated with the domain of this operator. More...
 
virtual const Epetra_MapOperatorRangeMap () const
 Returns the Epetra_Map object associated with the range of this operator (same as domain). More...
 
virtual const Epetra_BlockMapMap () const
 Implement the Epetra_SrcDistObjec::Map() function. More...
 
virtual const Epetra_MapRowMatrixRowMap () const
 Returns the Row Map object needed for implementing Epetra_RowMatrix. More...
 
virtual 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...
 
virtual const Epetra_CommComm () const
 Returns a pointer to the Epetra_Comm communicator associated with this matrix. More...
 
virtual void Print (std::ostream &os) const
 Print method. More...
 
virtual int SetUseTranspose (bool use_transpose)
 If set true, transpose of this operator will be applied. More...
 
virtual const char * Label () const
 Returns a character string describing the operator. More...
 
virtual int Apply (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_RowMatrix applied to a Epetra_MultiVector X in Y. More...
 
virtual int ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const
 Returns the result of a Epetra_RowMatrix inverse applied to an Epetra_MultiVector X in Y. More...
 
bool HasNormInf () const
 Returns true because this class can compute an Inf-norm. More...
 
virtual bool UseTranspose () const
 Returns the current UseTranspose setting. More...
 
virtual const Epetra_ImportImporter () const
 Returns the Epetra_Import object that contains the import operations for distributed operations, returns zero if none. More...
 
virtual const Epetra_ExportExporter () const
 Returns the Epetra_Export object that contains the export operations for distributed operations, returns zero if none. More...
 
- Public Member Functions inherited from Epetra_CompObject
Epetra_CompObjectoperator= (const Epetra_CompObject &src)
 
 Epetra_CompObject ()
 Basic Epetra_CompObject constuctor. More...
 
 Epetra_CompObject (const Epetra_CompObject &Source)
 Epetra_CompObject copy constructor. More...
 
virtual ~Epetra_CompObject ()
 Epetra_CompObject destructor. More...
 
void SetFlopCounter (const Epetra_Flops &FlopCounter_in)
 Set the internal Epetra_Flops() pointer. More...
 
void SetFlopCounter (const Epetra_CompObject &CompObject)
 Set the internal Epetra_Flops() pointer to the flop counter of another Epetra_CompObject. More...
 
void UnsetFlopCounter ()
 Set the internal Epetra_Flops() pointer to 0 (no flops counted). More...
 
Epetra_FlopsGetFlopCounter () const
 Get the pointer to the Epetra_Flops() object associated with this object, returns 0 if none. More...
 
void ResetFlops () const
 Resets the number of floating point operations to zero for this multi-vector. More...
 
double Flops () const
 Returns the number of floating point operations with this multi-vector. More...
 
void UpdateFlops (int Flops_in) const
 Increment Flop count for this object. More...
 
void UpdateFlops (long int Flops_in) const
 Increment Flop count for this object. More...
 
void UpdateFlops (long long Flops_in) const
 Increment Flop count for this object. More...
 
void UpdateFlops (double Flops_in) const
 Increment Flop count for this object. More...
 
void UpdateFlops (float Flops_in) const
 Increment Flop count for this object. More...
 
- Public Member Functions inherited from Epetra_Object
 Epetra_Object (int TracebackModeIn=-1, bool set_label=true)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const char *const Label, int TracebackModeIn=-1)
 Epetra_Object Constructor. More...
 
 Epetra_Object (const Epetra_Object &Object)
 Epetra_Object Copy Constructor. More...
 
virtual ~Epetra_Object ()
 Epetra_Object Destructor. More...
 
virtual void SetLabel (const char *const Label)
 Epetra_Object Label definition using char *. More...
 
virtual int ReportError (const std::string Message, int ErrorCode) const
 Error reporting method. More...
 
- Public Member Functions inherited from Epetra_RowMatrix
virtual ~Epetra_RowMatrix ()
 Destructor. More...
 
- Public Member Functions inherited from Epetra_Operator
virtual ~Epetra_Operator ()
 Destructor. More...
 
- Public Member Functions inherited from Epetra_SrcDistObject
virtual ~Epetra_SrcDistObject ()
 Epetra_SrcDistObject destructor. More...
 
- Static Public Member Functions inherited from Epetra_Object
static void SetTracebackMode (int TracebackModeValue)
 Set the value of the Epetra_Object error traceback report mode. More...
 
static int GetTracebackMode ()
 Get the value of the Epetra_Object error report mode. More...
 
static std::ostream & GetTracebackStream ()
 Get the output stream for error reporting. More...
 
- Static Public Attributes inherited from Epetra_Object
static int TracebackMode
 
- Protected Member Functions inherited from Epetra_BasicRowMatrix
void Setup ()
 
void UpdateImportVector (int NumVectors) const
 
void UpdateExportVector (int NumVectors) const
 
void SetImportExport ()
 
virtual void ComputeStructureConstants () const
 Update the constants associated with the structure of the matrix: Call only if structure changes from the initial RowMatrix. More...
 
virtual void ComputeNumericConstants () const
 Update the constants associated with the values of the matrix: Call only if values changes from the initial RowMatrix. More...
 
- Protected Member Functions inherited from Epetra_Object
std::string toString (const int &x) const
 
std::string toString (const long long &x) const
 
std::string toString (const double &x) const
 
- Protected Attributes inherited from Epetra_BasicRowMatrix
Epetra_CommComm_
 
Epetra_Map OperatorDomainMap_
 
Epetra_Map OperatorRangeMap_
 
Epetra_Map RowMatrixRowMap_
 
Epetra_Map RowMatrixColMap_
 
int NumMyNonzeros_
 
long long NumGlobalNonzeros_
 
int MaxNumEntries_
 
double NormInf_
 
double NormOne_
 
int NumMyRows_
 
int NumMyCols_
 
bool UseTranspose_
 
bool HasNormInf_
 
bool LowerTriangular_
 
bool UpperTriangular_
 
bool HaveStructureConstants_
 
bool HaveNumericConstants_
 
bool HaveMaps_
 
Epetra_MultiVectorImportVector_
 
Epetra_MultiVectorExportVector_
 
Epetra_ImportImporter_
 
Epetra_ExportExporter_
 
- Protected Attributes inherited from Epetra_CompObject
Epetra_FlopsFlopCounter_
 

Detailed Description

Epetra_VbrRowMatrix: A class for using an existing Epetra_VbrMatrix object as an Epetra_RowMatrix object.

The Epetra_VbrRowMatrix class takes an existing Epetra_VbrMatrix object and allows its use as an Epetra_RowMatrix without allocating additional storage. Although the Epetra_VbrMatrix itself inherits from Epetra_RowMatrix, a design flaw in the inheritance structure of Epetra prohibits the use of an Epetra_VbrMatrix object as an Epetra_RowMatrix in some important situations. Therefore we recommend the use of this class to wrap an Epetra_VbrMatrix object.

Warning
This class takes a pointer to an existing Epetra_VbrMatrix object. It is assumed that the user will pass in a pointer to a valid Epetra_VbrMatrix object, and will retain it throughout the life of the Epetra_VbrRowMatrix object.

Definition at line 68 of file Epetra_VbrRowMatrix.h.

Constructor & Destructor Documentation

Epetra_VbrRowMatrix::Epetra_VbrRowMatrix ( Epetra_VbrMatrix Matrix)
inline

Epetra_VbrRowMatrix constuctor.

Definition at line 81 of file Epetra_VbrRowMatrix.h.

virtual Epetra_VbrRowMatrix::~Epetra_VbrRowMatrix ( )
inlinevirtual

Epetra_VbrRowMatrix Destructor.

Definition at line 89 of file Epetra_VbrRowMatrix.h.

Member Function Documentation

int Epetra_VbrRowMatrix::UpdateMatrix ( Epetra_VbrMatrix Matrix)
inline

Update the matrix to which this object points.

Definition at line 99 of file Epetra_VbrRowMatrix.h.

int Epetra_VbrRowMatrix::ExtractMyRowCopy ( int  MyRow,
int  Length,
int &  NumEntries,
double *  Values,
int *  Indices 
) const
inlinevirtual

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

Parameters
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.
Returns
Integer error code, set to 0 if successful, set to -1 if MyRow not valid, -2 if Length is too short (NumEntries will have required length).

Implements Epetra_BasicRowMatrix.

Definition at line 121 of file Epetra_VbrRowMatrix.h.

int Epetra_VbrRowMatrix::ExtractMyEntryView ( int  CurEntry,
double *&  Value,
int &  RowIndex,
int &  ColIndex 
)
inlinevirtual

Returns a reference to the ith entry in the matrix, along with its row and column index.

Parameters
CurEntry(In) - Local entry to extract.
Value(Out) - Extracted reference to current values.
RowIndex(Out) - Row index for current entry.
ColIndex(Out) - Column index for current entry.
Returns
Integer error code, set to 0 if successful, set to -1 if CurEntry not valid.

Implements Epetra_BasicRowMatrix.

Definition at line 136 of file Epetra_VbrRowMatrix.h.

int Epetra_VbrRowMatrix::ExtractMyEntryView ( int  CurEntry,
double const *&  Value,
int &  RowIndex,
int &  ColIndex 
) const
inlinevirtual

Returns a const reference to the ith entry in the matrix, along with its row and column index.

Parameters
CurEntry(In) - Local entry to extract.
Value(Out) - Extracted reference to current values.
RowIndex(Out) - Row index for current entry.
ColIndex(Out) - Column index for current entry.
Returns
Integer error code, set to 0 if successful, set to -1 if CurEntry not valid.

Implements Epetra_BasicRowMatrix.

Definition at line 149 of file Epetra_VbrRowMatrix.h.

int Epetra_VbrRowMatrix::NumMyRowEntries ( int  MyRow,
int &  NumEntries 
) const
inlinevirtual

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
MyRow- (In) Local row.
NumEntries- (Out) Number of nonzero values.
Returns
Integer error code, set to 0 if successful, set to -1 if MyRow not valid.
Precondition
None.
Postcondition
Unchanged.

Implements Epetra_BasicRowMatrix.

Definition at line 163 of file Epetra_VbrRowMatrix.h.

int Epetra_VbrRowMatrix::RightScale ( const Epetra_Vector x)
inlinevirtual

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

Reimplemented from Epetra_BasicRowMatrix.

Definition at line 181 of file Epetra_VbrRowMatrix.h.

int Epetra_VbrRowMatrix::LeftScale ( const Epetra_Vector x)
inlinevirtual

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

Reimplemented from Epetra_BasicRowMatrix.

Definition at line 196 of file Epetra_VbrRowMatrix.h.

int Epetra_VbrRowMatrix::Multiply ( bool  TransA,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

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

Reimplemented from Epetra_BasicRowMatrix.

Definition at line 214 of file Epetra_VbrRowMatrix.h.

int Epetra_VbrRowMatrix::Solve ( bool  Upper,
bool  Trans,
bool  UnitDiagonal,
const Epetra_MultiVector X,
Epetra_MultiVector Y 
) const
inlinevirtual

Returns the result of a Epetra_VbrRowMatrix solve with a Epetra_MultiVector X in Y (not implemented).

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.

Reimplemented from Epetra_BasicRowMatrix.

Definition at line 234 of file Epetra_VbrRowMatrix.h.

Member Data Documentation

Epetra_VbrMatrix* Epetra_VbrRowMatrix::matrix_
private

Definition at line 245 of file Epetra_VbrRowMatrix.h.


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