10 #include "Teko_BlockedTpetraOperator.hpp"
11 #include "Teko_TpetraBlockedMappingStrategy.hpp"
12 #include "Teko_TpetraReorderedMappingStrategy.hpp"
14 #include "Teuchos_VerboseObject.hpp"
16 #include "Thyra_LinearOpBase.hpp"
17 #include "Thyra_TpetraLinearOp.hpp"
18 #include "Thyra_TpetraThyraWrappers.hpp"
19 #include "Thyra_DefaultProductMultiVector.hpp"
20 #include "Thyra_DefaultProductVectorSpace.hpp"
21 #include "Thyra_DefaultBlockedLinearOp.hpp"
23 #include "MatrixMarket_Tpetra.hpp"
28 namespace TpetraHelpers {
32 using Teuchos::rcp_dynamic_cast;
35 const std::vector<std::vector<GO> >& vars,
36 const Teuchos::RCP<
const Tpetra::Operator<ST, LO, GO, NT> >& content,
const std::string& label)
42 const std::vector<std::vector<GO> >& vars,
43 const Teuchos::RCP<
const Tpetra::Operator<ST, LO, GO, NT> >& content) {
44 fullContent_ = content;
45 blockedMapping_ = rcp(
new TpetraBlockedMappingStrategy(vars, fullContent_->getDomainMap(),
46 *fullContent_->getDomainMap()->getComm()));
47 SetMapStrategy(blockedMapping_);
50 BuildBlockedOperator();
53 void BlockedTpetraOperator::BuildBlockedOperator() {
54 TEUCHOS_ASSERT(blockedMapping_ != Teuchos::null);
57 const RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> > crsContent =
58 rcp_dynamic_cast<
const Tpetra::CrsMatrix<ST, LO, GO, NT> >(fullContent_);
61 if (blockedOperator_ == Teuchos::null) {
62 blockedOperator_ = blockedMapping_->buildBlockedThyraOp(crsContent, label_);
64 const RCP<Thyra::BlockedLinearOpBase<ST> > blkOp =
65 rcp_dynamic_cast<Thyra::BlockedLinearOpBase<ST> >(blockedOperator_,
true);
66 blockedMapping_->rebuildBlockedThyraOp(crsContent, blkOp);
70 SetOperator(blockedOperator_,
false);
73 if (reorderManager_ != Teuchos::null)
Reorder(*reorderManager_);
76 const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > BlockedTpetraOperator::GetBlock(
78 const RCP<const Thyra::BlockedLinearOpBase<ST> > blkOp =
79 Teuchos::rcp_dynamic_cast<
const Thyra::BlockedLinearOpBase<ST> >(
getThyraOp());
81 RCP<const Thyra::TpetraLinearOp<ST, LO, GO, NT> > tOp =
82 rcp_dynamic_cast<
const Thyra::TpetraLinearOp<ST, LO, GO, NT> >(blkOp->getBlock(i, j),
true);
83 return tOp->getConstTpetraOperator();
93 RCP<const MappingStrategy> reorderMapping =
94 rcp(
new TpetraReorderedMappingStrategy(*reorderManager_, blockedMapping_));
95 RCP<const Thyra::BlockedLinearOpBase<ST> > blockOp =
96 rcp_dynamic_cast<
const Thyra::BlockedLinearOpBase<ST> >(blockedOperator_);
98 RCP<const Thyra::LinearOpBase<ST> > A = buildReorderedLinearOp(*reorderManager_, blockOp);
101 SetMapStrategy(reorderMapping);
102 SetOperator(A,
false);
107 SetMapStrategy(blockedMapping_);
108 SetOperator(blockedOperator_,
false);
109 reorderManager_ = Teuchos::null;
115 RCP<Thyra::PhysicallyBlockedLinearOpBase<ST> > blockOp =
116 rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<ST> >(blockedOperator_);
119 int rows = Teko::blockRowCount(blockOp);
121 for (
int i = 0; i < rows; i++) {
122 for (
int j = 0; j < rows; j++) {
124 RCP<const Thyra::TpetraLinearOp<ST, LO, GO, NT> > tOp =
125 rcp_dynamic_cast<
const Thyra::TpetraLinearOp<ST, LO, GO, NT> >(blockOp->getBlock(i, j));
126 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> > mat =
127 Teuchos::rcp_dynamic_cast<
const Tpetra::CrsMatrix<ST, LO, GO, NT> >(
128 tOp->getConstTpetraOperator());
131 Tpetra::MatrixMarket::Writer<Tpetra::CrsMatrix<ST, LO, GO, NT> >::writeSparseFile(
132 formatBlockName(prefix, i, j, rows).c_str(), mat);
138 Tpetra::Vector<ST, LO, GO, NT> xf(getRangeMap());
139 Tpetra::Vector<ST, LO, GO, NT> xs(getRangeMap());
140 Tpetra::Vector<ST, LO, GO, NT> y(getDomainMap());
144 ST diffNorm = 0.0, trueNorm = 0.0;
145 for (
int i = 0; i < count; i++) {
152 fullContent_->apply(y, xf);
155 xs.update(-1.0, xf, 1.0);
156 diffNorm = xs.norm2();
157 trueNorm = xf.norm2();
160 result &= (diffNorm / trueNorm < tol);
void Reorder(const BlockReorderManager &brm)
Class that describes how a flat blocked operator should be reordered.
void RemoveReording()
Remove any reordering on this object.
const RCP< const Thyra::LinearOpBase< ST > > getThyraOp() const
Return the thyra operator associated with this wrapper.
virtual void WriteBlocks(const std::string &prefix) const
BlockedTpetraOperator(const std::vector< std::vector< GO > > &vars, const Teuchos::RCP< const Tpetra::Operator< ST, LO, GO, NT > > &content, const std::string &label="<ANYM>")
bool testAgainstFullOperator(int count, ST tol) const
Helps perform sanity checks.
virtual void SetContent(const std::vector< std::vector< GO > > &vars, const Teuchos::RCP< const Tpetra::Operator< ST, LO, GO, NT > > &content)
Implements the Epetra_Operator interface with a Thyra LinearOperator. This enables the use of absrtac...