NOX
Development
|
An Epetra row matrix for implementing the operator . More...
#include <LOCA_Epetra_LowRankUpdateRowMatrix.H>
Public Member Functions | |
LowRankUpdateRowMatrix (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< Epetra_RowMatrix > &jacRowMatrix, const Teuchos::RCP< Epetra_MultiVector > &U_multiVec, const Teuchos::RCP< Epetra_MultiVector > &V_multiVec, bool setup_for_solve, bool include_UV_terms) | |
Constructor. More... | |
virtual | ~LowRankUpdateRowMatrix () |
Destructor. | |
virtual const Epetra_BlockMap & | Map () const |
Returns a reference to the Epetra_BlockMap for this object. | |
virtual int | NumMyRowEntries (int MyRow, int &NumEntries) const |
Returns the number of nonzero entries in MyRow. | |
virtual int | MaxNumEntries () const |
Returns the maximum of NumMyRowEntries() over all rows. | |
virtual 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. | |
virtual int | ExtractDiagonalCopy (Epetra_Vector &Diagonal) const |
Returns a copy of the main diagonal in a user-provided vector. | |
virtual int | Multiply (bool TransA, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Returns the result of a Epetra_RowMatrix multiplied by a Epetra_MultiVector X in Y. | |
virtual int | Solve (bool Upper, bool Trans, bool UnitDiagonal, const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
Returns result of a local-only solve using a triangular Epetra_RowMatrix with Epetra_MultiVectors X and Y. | |
virtual int | InvRowSums (Epetra_Vector &x) const |
Computes the sum of absolute values of the rows of the Epetra_RowMatrix, results returned in x. | |
virtual int | LeftScale (const Epetra_Vector &x) |
Scales the Epetra_RowMatrix on the left with a Epetra_Vector x. | |
virtual int | InvColSums (Epetra_Vector &x) const |
Computes the sum of absolute values of the columns of the Epetra_RowMatrix, results returned in x. | |
virtual int | RightScale (const Epetra_Vector &x) |
Scales the Epetra_RowMatrix on the right with a Epetra_Vector x. | |
virtual bool | Filled () const |
If FillComplete() has been called, this query returns true, otherwise it returns false. | |
virtual double | NormInf () const |
Returns the infinity norm of the global matrix. | |
virtual double | NormOne () const |
Returns the one norm of the global matrix. | |
virtual int | NumGlobalNonzeros () const |
Returns the number of nonzero entries in the global matrix. | |
virtual long long | NumGlobalNonzeros64 () const |
virtual int | NumGlobalRows () const |
Returns the number of global matrix rows. | |
virtual long long | NumGlobalRows64 () const |
virtual int | NumGlobalCols () const |
Returns the number of global matrix columns. | |
virtual long long | NumGlobalCols64 () const |
virtual int | NumGlobalDiagonals () const |
Returns the number of global nonzero diagonal entries, based on global row/column index comparisons. | |
virtual long long | NumGlobalDiagonals64 () const |
virtual int | NumMyNonzeros () const |
Returns the number of nonzero entries in the calling processor's portion of the matrix. | |
virtual int | NumMyRows () const |
Returns the number of matrix rows owned by the calling processor. | |
virtual int | NumMyCols () const |
Returns the number of matrix columns owned by the calling processor. | |
virtual int | NumMyDiagonals () const |
Returns the number of local nonzero diagonal entries, based on global row/column index comparisons. | |
virtual bool | LowerTriangular () const |
If matrix is lower triangular in local index space, this query returns true, otherwise it returns false. | |
virtual bool | UpperTriangular () const |
If matrix is upper triangular in local index space, this query returns true, otherwise it returns false. | |
virtual const Epetra_Map & | RowMatrixRowMap () const |
Returns the Epetra_Map object associated with the rows of this matrix. | |
virtual const Epetra_Map & | RowMatrixColMap () const |
Returns the Epetra_Map object associated with the columns of this matrix. | |
virtual const Epetra_Import * | RowMatrixImporter () const |
Returns the Epetra_Import object that contains the import operations for distributed operations. | |
Public Member Functions inherited from LOCA::Epetra::LowRankUpdateOp | |
LowRankUpdateOp (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< Epetra_Operator > &jacOperator, const Teuchos::RCP< const Epetra_MultiVector > &U_multiVec, const Teuchos::RCP< const Epetra_MultiVector > &V_multiVec, bool setup_for_solve) | |
Constructor. More... | |
virtual | ~LowRankUpdateOp () |
Destructor. | |
virtual int | SetUseTranspose (bool UseTranspose) |
Set to true if the transpose of the operator is requested. | |
virtual int | Apply (const Epetra_MultiVector &Input, Epetra_MultiVector &Result) const |
Returns the result of a Epetra_Operator applied to a Epetra_MultiVector Input in Result as described above. | |
virtual int | ApplyInverse (const Epetra_MultiVector &X, Epetra_MultiVector &Y) const |
This method does nothing. | |
virtual const char * | Label () const |
Returns a character std::string describing the operator. | |
virtual bool | UseTranspose () const |
Returns the current UseTranspose setting. Always returns false. | |
virtual bool | HasNormInf () const |
Returns true if the this object can provide an approximate Inf-norm, false otherwise. | |
virtual const Epetra_Comm & | Comm () const |
Returns a reference to the Epetra_Comm communicator associated with this operator. | |
virtual const Epetra_Map & | OperatorDomainMap () const |
Returns the Epetra_Map object associated with the domain of this matrix operator. | |
virtual const Epetra_Map & | OperatorRangeMap () const |
Returns the Epetra_Map object associated with the range of this matrix operator. | |
Protected Member Functions | |
double | computeUV (int MyRow, int MyCol) const |
Compute MyRow , MyCol entry of . | |
Protected Attributes | |
Teuchos::RCP< Epetra_RowMatrix > | J_rowMatrix |
Stores row matrix representing J. | |
Teuchos::RCP< Epetra_MultiVector > | nonconst_U |
Stores pointer to non-const U. | |
Teuchos::RCP< Epetra_MultiVector > | nonconst_V |
Stores pointer to non-const V. | |
bool | includeUV |
Flag indicating whether to include U*V^T terms. | |
int | m |
Number of columns in U and V. | |
const Epetra_BlockMap & | U_map |
Map for U. | |
const Epetra_BlockMap & | V_map |
Map for V. | |
const Epetra_BlockMap & | row_map |
Row map for J. | |
Protected Attributes inherited from LOCA::Epetra::LowRankUpdateOp | |
Teuchos::RCP< LOCA::GlobalData > | globalData |
Global data object. | |
std::string | label |
Label for operator. | |
Epetra_LocalMap | localMap |
Local map for generating temporary matrices. | |
Teuchos::RCP< Epetra_Operator > | J |
Stores operator representing J. | |
Teuchos::RCP< const Epetra_MultiVector > | U |
Stores multivector representing U. | |
Teuchos::RCP< const Epetra_MultiVector > | V |
Stores multivector representing V. | |
bool | useTranspose |
Flag indicating whether to use the transpose. | |
Teuchos::RCP< Epetra_MultiVector > | tmpMat |
Temporary matrix. | |
Teuchos::RCP< Epetra_MultiVector > | JinvU |
Stores J^{-1}*U. | |
Teuchos::RCP< Epetra_MultiVector > | lu |
Stores LU factorization of I + V^T*J^{-1}*U. | |
std::vector< int > | ipiv |
Stores pivots for LU factorization. | |
Teuchos::LAPACK< int, double > | lapack |
Lapack wrappers. | |
An Epetra row matrix for implementing the operator .
This class implements the Epetra_RowMatrix interface for where is an Epetra_RowMatrix and and are Epetra_MultiVectors. It is derived from LOCA::Epetra::LowRankUpdateOp to implement the Epetra_Operator interface. The interface here implements the Epetra_RowMatrix interface when the matrix is itself a row matrix. This allows preconditioners to be computed and scaling in linear systems to be performed when using this operator. The implementation here merely adds the corresponding entries for to the rows of . Note however this is only an approximation to the true matrix .
This class assumes and have the same distribution as the rows of .
LOCA::Epetra::LowRankUpdateRowMatrix::LowRankUpdateRowMatrix | ( | const Teuchos::RCP< LOCA::GlobalData > & | global_data, |
const Teuchos::RCP< Epetra_RowMatrix > & | jacRowMatrix, | ||
const Teuchos::RCP< Epetra_MultiVector > & | U_multiVec, | ||
const Teuchos::RCP< Epetra_MultiVector > & | V_multiVec, | ||
bool | setup_for_solve, | ||
bool | include_UV_terms | ||
) |
Constructor.
global_data | [in] The global data object |
jacRowMatrix | [in] Jacobian operator J as a row matrix |
U_multiVec | [in] Multivector representing U |
V_multiVec | [in] Multivector representing V |
setup_for_solve | [in] Setup data structures for ApplyInverse() |
include_UV_terms | [in] Include terms in RowMatrix routines ExtractRowCopy(), ExtactDiagonalCopy(), InvRowSums(), InvColSums(), NormInf() and NormOne(). |