NOX
Development
|
A Tpetra row matrix for implementing the operator . More...
#include <LOCA_Tpetra_LowRankUpdateRowMatrix.hpp>
Public Member Functions | |
LowRankUpdateRowMatrix (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< NOX::TRowMatrix > &jacRowMatrix, const Teuchos::RCP< NOX::TMultiVector > &U_multiVec, const Teuchos::RCP< NOX::TMultiVector > &V_multiVec, bool setup_for_solve, bool include_UV_terms) | |
Constructor. More... | |
virtual local_ordinal_type | getBlockSize () const override |
virtual Teuchos::RCP< const Teuchos::Comm< int > > | getComm () const override |
virtual Teuchos::RCP< const NOX::TMap > | getRowMap () const override |
virtual Teuchos::RCP< const NOX::TMap > | getColMap () const override |
virtual Teuchos::RCP< const NOX::TRowGraph > | getGraph () const override |
virtual ::Tpetra::global_size_t | getGlobalNumRows () const override |
virtual ::Tpetra::global_size_t | getGlobalNumCols () const override |
virtual size_t | getLocalNumRows () const override |
virtual size_t | getLocalNumCols () const override |
virtual NOX::GlobalOrdinal | getIndexBase () const override |
virtual ::Tpetra::global_size_t | getGlobalNumEntries () const override |
virtual size_t | getLocalNumEntries () const override |
virtual size_t | getNumEntriesInGlobalRow (NOX::GlobalOrdinal globalRow) const override |
virtual size_t | getNumEntriesInLocalRow (NOX::LocalOrdinal localRow) const override |
virtual size_t | getGlobalMaxNumRowEntries () const override |
virtual size_t | getLocalMaxNumRowEntries () const override |
virtual bool | hasColMap () const override |
virtual bool | isLocallyIndexed () const override |
virtual bool | isGloballyIndexed () const override |
virtual bool | isFillComplete () const override |
virtual bool | supportsRowViews () const override |
virtual void | getGlobalRowCopy (NOX::GlobalOrdinal GlobalRow, NOX::TRowMatrix::nonconst_global_inds_host_view_type &Indices, NOX::TRowMatrix::nonconst_values_host_view_type &Values, size_t &NumEntries) const override |
virtual void | getLocalRowCopy (NOX::LocalOrdinal LocalRow, NOX::TRowMatrix::nonconst_local_inds_host_view_type &Indices, NOX::TRowMatrix::nonconst_values_host_view_type &Values, size_t &NumEntries) const override |
virtual void | getGlobalRowView (NOX::GlobalOrdinal GlobalRow, NOX::TRowMatrix::global_inds_host_view_type &Indices, NOX::TRowMatrix::values_host_view_type &Values) const override |
virtual void | getLocalRowView (NOX::LocalOrdinal LocalRow, NOX::TRowMatrix::local_inds_host_view_type &Indices, NOX::TRowMatrix::values_host_view_type &Values) const override |
virtual void | getLocalDiagCopy (NOX::TVector &diag) const override |
virtual void | leftScale (const NOX::TVector &x) override |
virtual void | rightScale (const NOX::TVector &x) override |
virtual mag_type | getFrobeniusNorm () const override |
virtual Teuchos::RCP< const NOX::TMap > | getDomainMap () const override |
virtual Teuchos::RCP< const NOX::TMap > | getRangeMap () const override |
virtual void | apply (const NOX::TMultiVector &X, NOX::TMultiVector &Y, Teuchos::ETransp mode=Teuchos::NO_TRANS, NOX::Scalar alpha=Teuchos::ScalarTraits< NOX::Scalar >::one(), NOX::Scalar beta=Teuchos::ScalarTraits< NOX::Scalar >::zero()) const override |
Protected Member Functions | |
template<typename ViewType > | |
KOKKOS_INLINE_FUNCTION NOX::Scalar | computeUV (int MyRow, int MyCol) const |
Compute MyRow , MyCol entry of . Views are local Kokkos view types. | |
Protected Attributes | |
Teuchos::RCP< NOX::TRowMatrix > | J_rowMatrix |
Stores row matrix representing J. | |
Teuchos::RCP< NOX::TMultiVector > | nonconst_U |
Stores pointer to non-const U. | |
Teuchos::RCP< NOX::TMultiVector > | 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 NOX::TMap & | U_map |
Map for U. | |
const NOX::TMap & | V_map |
Map for V. | |
const NOX::TMap & | row_map |
Row map for J. | |
Teuchos::RCP< NOX::TMultiVector > | tmpMat |
Temporary workspace. | |
Teuchos::RCP< NOX::TMap > | local_map |
Locally replicated map. | |
A Tpetra row matrix for implementing the operator .
This class implements the Tpetra::RowMatrix interface for where is an Tpetra::RowMatrix and and are Tpetra::MultiVectors. It is derived from LOCA::Tpetra::LowRankUpdateOp to implement the Tpetra::Operator interface. The interface here implements the Tpetra::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 (which is dense).
This class assumes and have the same distribution as the rows of .
LOCA::Tpetra::LowRankUpdateRowMatrix::LowRankUpdateRowMatrix | ( | const Teuchos::RCP< LOCA::GlobalData > & | global_data, |
const Teuchos::RCP< NOX::TRowMatrix > & | jacRowMatrix, | ||
const Teuchos::RCP< NOX::TMultiVector > & | U_multiVec, | ||
const Teuchos::RCP< NOX::TMultiVector > & | 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(). |
References local_map, and Teuchos::rcp().