10 #include "Teko_BlockedEpetraOperator.hpp"
11 #include "Teko_BlockedMappingStrategy.hpp"
12 #include "Teko_ReorderedMappingStrategy.hpp"
14 #include "Teuchos_VerboseObject.hpp"
16 #include "Thyra_LinearOpBase.hpp"
17 #include "Thyra_EpetraLinearOp.hpp"
18 #include "Thyra_EpetraThyraWrappers.hpp"
19 #include "Thyra_DefaultProductMultiVector.hpp"
20 #include "Thyra_DefaultProductVectorSpace.hpp"
21 #include "Thyra_DefaultBlockedLinearOp.hpp"
22 #include "Thyra_get_Epetra_Operator.hpp"
24 #include "EpetraExt_MultiVectorOut.h"
25 #include "EpetraExt_RowMatrixOut.h"
34 using Teuchos::rcp_dynamic_cast;
37 const Teuchos::RCP<const Epetra_Operator>& content,
38 const std::string& label)
44 const Teuchos::RCP<const Epetra_Operator>& content) {
45 fullContent_ = content;
46 blockedMapping_ = rcp(
new BlockedMappingStrategy(
47 vars, Teuchos::rcpFromRef(fullContent_->OperatorDomainMap()), fullContent_->Comm()));
48 SetMapStrategy(blockedMapping_);
51 BuildBlockedOperator();
54 void BlockedEpetraOperator::BuildBlockedOperator() {
55 TEUCHOS_ASSERT(blockedMapping_ != Teuchos::null);
58 const RCP<const Epetra_CrsMatrix> crsContent =
59 rcp_dynamic_cast<
const Epetra_CrsMatrix>(fullContent_);
62 if (blockedOperator_ == Teuchos::null) {
63 blockedOperator_ = blockedMapping_->buildBlockedThyraOp(crsContent, label_);
65 const RCP<Thyra::BlockedLinearOpBase<double> > blkOp =
66 rcp_dynamic_cast<Thyra::BlockedLinearOpBase<double> >(blockedOperator_,
true);
67 blockedMapping_->rebuildBlockedThyraOp(crsContent, blkOp);
71 SetOperator(blockedOperator_,
false);
74 if (reorderManager_ != Teuchos::null)
Reorder(*reorderManager_);
77 const Teuchos::RCP<const Epetra_Operator> BlockedEpetraOperator::GetBlock(
int i,
int j)
const {
78 const RCP<const Thyra::BlockedLinearOpBase<double> > blkOp =
79 Teuchos::rcp_dynamic_cast<
const Thyra::BlockedLinearOpBase<double> >(
getThyraOp());
81 return Thyra::get_Epetra_Operator(*blkOp->getBlock(i, j));
91 RCP<const MappingStrategy> reorderMapping =
92 rcp(
new ReorderedMappingStrategy(*reorderManager_, blockedMapping_));
93 RCP<const Thyra::BlockedLinearOpBase<double> > blockOp =
94 rcp_dynamic_cast<
const Thyra::BlockedLinearOpBase<double> >(blockedOperator_);
96 RCP<const Thyra::LinearOpBase<double> > A = buildReorderedLinearOp(*reorderManager_, blockOp);
99 SetMapStrategy(reorderMapping);
100 SetOperator(A,
false);
105 SetMapStrategy(blockedMapping_);
106 SetOperator(blockedOperator_,
false);
107 reorderManager_ = Teuchos::null;
113 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blockOp =
114 rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<double> >(blockedOperator_);
117 int rows = Teko::blockRowCount(blockOp);
119 for (
int i = 0; i < rows; i++) {
120 for (
int j = 0; j < rows; j++) {
122 RCP<const Epetra_RowMatrix> mat = Teuchos::rcp_dynamic_cast<
const Epetra_RowMatrix>(
123 Thyra::get_Epetra_Operator(*blockOp->getBlock(i, j)));
126 EpetraExt::RowMatrixToMatrixMarketFile(formatBlockName(prefix, i, j, rows).c_str(), *mat);
132 Epetra_Vector xf(OperatorRangeMap());
133 Epetra_Vector xs(OperatorRangeMap());
134 Epetra_Vector y(OperatorDomainMap());
138 double diffNorm = 0.0, trueNorm = 0.0;
139 for (
int i = 0; i < count; i++) {
146 fullContent_->Apply(y, xf);
149 xs.Update(-1.0, xf, 1.0);
154 result &= (diffNorm / trueNorm < tol);
bool testAgainstFullOperator(int count, double tol) const
Helps perform sanity checks.
Class that describes how a flat blocked operator should be reordered.
void Reorder(const BlockReorderManager &brm)
Implements the Epetra_Operator interface with a Thyra LinearOperator. This enables the use of absrtac...
virtual void SetContent(const std::vector< std::vector< int > > &vars, const Teuchos::RCP< const Epetra_Operator > &content)
void RemoveReording()
Remove any reordering on this object.
BlockedEpetraOperator(const std::vector< std::vector< int > > &vars, const Teuchos::RCP< const Epetra_Operator > &content, const std::string &label="<ANYM>")
const RCP< const Thyra::LinearOpBase< double > > getThyraOp() const
Return the thyra operator associated with this wrapper.
virtual void WriteBlocks(const std::string &prefix) const