Epetra_MsrMatrix: A class for constructing and using real-valued double-precision sparse compressed row matrices. More...
#include <Epetra_MsrMatrix.h>
Public Member Functions | |
Constructors/Destructor | |
Epetra_MsrMatrix (int *proc_config, AZ_MATRIX *Amat) | |
Epetra_MsrMatrix constuctor using existing Aztec DMSR matrix. More... | |
virtual | ~Epetra_MsrMatrix () |
Epetra_MsrMatrix Destructor. | |
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_MsrMatrix 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_MsrMatrix 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_MsrMatrix, results returned in x. More... | |
int | LeftScale (const Epetra_Vector &x) |
Scales the Epetra_MsrMatrix 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_MsrMatrix, results returned in x. More... | |
int | RightScale (const Epetra_Vector &x) |
Scales the Epetra_MsrMatrix 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. | |
bool | LowerTriangular () const |
If matrix is lower triangular, this query returns true, otherwise it returns false. | |
bool | UpperTriangular () const |
If matrix is upper triangular, this query returns true, otherwise it returns false. | |
Attribute access functions | |
AZ_MATRIX * | Amat () const |
Returns a pointer to the Aztec Msr matrix used to create this object. | |
double | NormInf () const |
Returns the infinity norm of the global matrix. | |
double | NormOne () const |
Returns the one norm of the global matrix. | |
int | NumGlobalNonzeros () const |
Returns the number of nonzero entries in the global matrix. | |
int | NumGlobalRows () const |
Returns the number of global matrix rows. | |
int | NumGlobalCols () const |
Returns the number of global matrix columns. | |
int | NumGlobalDiagonals () const |
Returns the number of global nonzero diagonal entries. | |
long long | NumGlobalNonzeros64 () const |
Returns the number of nonzero entries in the global matrix. | |
long long | NumGlobalRows64 () const |
Returns the number of global matrix rows. | |
long long | NumGlobalCols64 () const |
Returns the number of global matrix columns. | |
long long | NumGlobalDiagonals64 () const |
Returns the number of global nonzero diagonal entries. | |
int | NumMyNonzeros () const |
Returns the number of nonzero entries in the calling processor's portion of the matrix. | |
int | NumMyRows () const |
Returns the number of matrix rows owned by the calling processor. | |
int | NumMyCols () const |
Returns the number of matrix columns owned by the calling processor. | |
int | NumMyDiagonals () const |
Returns the number of local nonzero diagonal entries. | |
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 (same as domain). | |
const Epetra_BlockMap & | Map () const |
Implement the Epetra_SrcDistObjec::Map() function. | |
const Epetra_Map & | RowMatrixRowMap () const |
Returns the Row Map object needed for implementing Epetra_RowMatrix. | |
const Epetra_Map & | RowMatrixColMap () const |
Returns the Column Map object needed for implementing Epetra_RowMatrix. | |
virtual const Epetra_Import * | RowMatrixImporter () const |
Returns the Epetra_Import object that contains the import operations for distributed operations. | |
const Epetra_Comm & | Comm () const |
Returns a pointer to the Epetra_Comm communicator associated with this matrix. | |
I/O Methods | |
virtual void | Print (std::ostream &os) const |
Print method. | |
Additional methods required to support the Epetra_Operator interface | |
const char * | Label () const |
Returns a character std::string describing the operator. | |
int | SetUseTranspose (bool use_transpose) |
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. | |
virtual bool | UseTranspose () const |
Returns the current UseTranspose setting. | |
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. | |
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 int | ReportError (const std::string Message, int ErrorCode) const |
Error reporting method. | |
virtual void | SetLabel (const char *const Label) |
Epetra_Object Label definition using char *. More... | |
Public Member Functions inherited from Epetra_CompObject | |
Epetra_CompObject & | operator= (const Epetra_CompObject &src) |
Epetra_CompObject () | |
Basic Epetra_CompObject constuctor. | |
Epetra_CompObject (const Epetra_CompObject &Source) | |
Epetra_CompObject copy constructor. | |
virtual | ~Epetra_CompObject () |
Epetra_CompObject destructor. | |
void | SetFlopCounter (const Epetra_Flops &FlopCounter_in) |
Set the internal Epetra_Flops() pointer. | |
void | SetFlopCounter (const Epetra_CompObject &CompObject) |
Set the internal Epetra_Flops() pointer to the flop counter of another Epetra_CompObject. | |
void | UnsetFlopCounter () |
Set the internal Epetra_Flops() pointer to 0 (no flops counted). | |
Epetra_Flops * | GetFlopCounter () const |
Get the pointer to the Epetra_Flops() object associated with this object, returns 0 if none. | |
void | ResetFlops () const |
Resets the number of floating point operations to zero for this multi-vector. | |
double | Flops () const |
Returns the number of floating point operations with this multi-vector. | |
void | UpdateFlops (int Flops_in) const |
Increment Flop count for this object. | |
void | UpdateFlops (long int Flops_in) const |
Increment Flop count for this object. | |
void | UpdateFlops (long long Flops_in) const |
Increment Flop count for this object. | |
void | UpdateFlops (double Flops_in) const |
Increment Flop count for this object. | |
void | UpdateFlops (float Flops_in) const |
Increment Flop count for this object. | |
Public Member Functions inherited from Epetra_RowMatrix | |
virtual | ~Epetra_RowMatrix () |
Destructor. | |
Public Member Functions inherited from Epetra_Operator | |
virtual | ~Epetra_Operator () |
Destructor. | |
Public Member Functions inherited from Epetra_SrcDistObject | |
virtual | ~Epetra_SrcDistObject () |
Epetra_SrcDistObject destructor. | |
Additional Inherited Members | |
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. | |
static std::ostream & | GetTracebackStream () |
Get the output stream for error reporting. | |
Static Public Attributes inherited from Epetra_Object | |
static int | TracebackMode |
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_CompObject | |
Epetra_Flops * | FlopCounter_ |
Epetra_MsrMatrix: A class for constructing and using real-valued double-precision sparse compressed row matrices.
The Epetra_MsrMatrix provides basic support for existing Aztec users who have an investment in the Aztec DMSR matrix format. A user may pass an existing Aztec DMSR matrix to the constructor for this class. The data from the DMSR matrix will not be copied. Thus, any changes the user makes to the DMSR matrix data will be reflected in the associated Epetra_MsrMatrix object.
Epetra_MsrMatrix::Epetra_MsrMatrix | ( | int * | proc_config, |
AZ_MATRIX * | Amat | ||
) |
Epetra_MsrMatrix constuctor using existing Aztec DMSR matrix.
Creates a Epetra_MsrMatrix object by encapsulating an existing Aztec DMSR matrix. The Aztec matrix must come in as an AZ_MATRIX pointer, and AZ_transform must have called. Also, the AZ_matrix_type must be AZ_MSR_MATRIX. (If the matrix is stored in Amat, this information is contained in Amat->data_org[AZ_matrix_type].)
In | Amat - A completely constructed Aztec DMSR matrix. |
In | proc_config - An Aztec array containing information about the parallel machine. |
|
inlinevirtual |
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector X in Y.
In | X - A Epetra_MultiVector of dimension NumVectors to multiply with matrix. |
Out | Y -A Epetra_MultiVector of dimension NumVectors containing result. |
Implements Epetra_Operator.
References Multiply(), and UseTranspose().
|
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_MsrMatrix::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.
In | X - A Epetra_MultiVector of dimension NumVectors to solve for. |
Out | Y -A Epetra_MultiVector of dimension NumVectors containing result. |
Implements Epetra_Operator.
|
virtual |
Returns a copy of the main diagonal in a user-provided vector.
Out | Diagonal - Extracted main diagonal. |
Implements Epetra_RowMatrix.
|
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 global column indices for the corresponding values. |
Implements Epetra_RowMatrix.
|
virtual |
Computes the sum of absolute values of the columns of the Epetra_MsrMatrix, 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. |
Implements Epetra_RowMatrix.
|
virtual |
Computes the sum of absolute values of the rows of the Epetra_MsrMatrix, 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. |
Implements Epetra_RowMatrix.
|
virtual |
Scales the Epetra_MsrMatrix 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. |
Implements Epetra_RowMatrix.
|
virtual |
Returns the result of a Epetra_MsrMatrix 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. |
Implements Epetra_RowMatrix.
Referenced by Apply().
|
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.
In | MyRow - Local row. |
Out | NumEntries - Number of nonzero values. |
Implements Epetra_RowMatrix.
|
virtual |
Scales the Epetra_MsrMatrix 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. |
Implements Epetra_RowMatrix.
|
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.
In | use_transpose -If true, multiply by the transpose of operator, otherwise just use operator. |
Implements Epetra_Operator.
|
virtual |
Returns the result of a Epetra_MsrMatrix multiplied by a Epetra_MultiVector X in Y.
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. |
Implements Epetra_RowMatrix.