Teko Version of the Day
Loading...
Searching...
No Matches
Teko_BlockedEpetraOperator.cpp
1// @HEADER
2// *****************************************************************************
3// Teko: A package for block and physics based preconditioning
4//
5// Copyright 2010 NTESS and the Teko contributors.
6// SPDX-License-Identifier: BSD-3-Clause
7// *****************************************************************************
8// @HEADER
9
10#include "Teko_BlockedEpetraOperator.hpp"
11#include "Teko_BlockedMappingStrategy.hpp"
12#include "Teko_ReorderedMappingStrategy.hpp"
13
14#include "Teuchos_VerboseObject.hpp"
15
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"
23
24#include "EpetraExt_MultiVectorOut.h"
25#include "EpetraExt_RowMatrixOut.h"
26
27#include "Teko_Utilities.hpp"
28
29namespace Teko {
30namespace Epetra {
31
32using Teuchos::RCP;
33using Teuchos::rcp;
34using Teuchos::rcp_dynamic_cast;
35
36BlockedEpetraOperator::BlockedEpetraOperator(const std::vector<std::vector<int> >& vars,
37 const Teuchos::RCP<const Epetra_Operator>& content,
38 const std::string& label)
39 : Teko::Epetra::EpetraOperatorWrapper(), label_(label) {
40 SetContent(vars, content);
41}
42
43void BlockedEpetraOperator::SetContent(const std::vector<std::vector<int> >& vars,
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_);
49
50 // build thyra operator
51 BuildBlockedOperator();
52}
53
54void BlockedEpetraOperator::BuildBlockedOperator() {
55 TEUCHOS_ASSERT(blockedMapping_ != Teuchos::null);
56
57 // get a CRS matrix
58 const RCP<const Epetra_CrsMatrix> crsContent =
59 rcp_dynamic_cast<const Epetra_CrsMatrix>(fullContent_);
60
61 // ask the strategy to build the Thyra operator for you
62 if (blockedOperator_ == Teuchos::null) {
63 blockedOperator_ = blockedMapping_->buildBlockedThyraOp(crsContent, label_);
64 } else {
65 const RCP<Thyra::BlockedLinearOpBase<double> > blkOp =
66 rcp_dynamic_cast<Thyra::BlockedLinearOpBase<double> >(blockedOperator_, true);
67 blockedMapping_->rebuildBlockedThyraOp(crsContent, blkOp);
68 }
69
70 // set whatever is returned
71 SetOperator(blockedOperator_, false);
72
73 // reorder if neccessary
74 if (reorderManager_ != Teuchos::null) Reorder(*reorderManager_);
75}
76
77const 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());
80
81 return Thyra::get_Epetra_Operator(*blkOp->getBlock(i, j));
82}
83
88 reorderManager_ = rcp(new BlockReorderManager(brm));
89
90 // build reordered objects
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_);
95
96 RCP<const Thyra::LinearOpBase<double> > A = buildReorderedLinearOp(*reorderManager_, blockOp);
97
98 // set them as working values
99 SetMapStrategy(reorderMapping);
100 SetOperator(A, false);
101}
102
105 SetMapStrategy(blockedMapping_);
106 SetOperator(blockedOperator_, false);
107 reorderManager_ = Teuchos::null;
108}
109
112void BlockedEpetraOperator::WriteBlocks(const std::string& prefix) const {
113 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blockOp =
114 rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<double> >(blockedOperator_);
115
116 // get size of blocked block operator
117 int rows = Teko::blockRowCount(blockOp);
118
119 for (int i = 0; i < rows; i++) {
120 for (int j = 0; j < rows; j++) {
121 // get the row matrix object
122 RCP<const Epetra_RowMatrix> mat = Teuchos::rcp_dynamic_cast<const Epetra_RowMatrix>(
123 Thyra::get_Epetra_Operator(*blockOp->getBlock(i, j)));
124
125 // write to file
126 EpetraExt::RowMatrixToMatrixMarketFile(formatBlockName(prefix, i, j, rows).c_str(), *mat);
127 }
128 }
129}
130
131bool BlockedEpetraOperator::testAgainstFullOperator(int count, double tol) const {
132 Epetra_Vector xf(OperatorRangeMap());
133 Epetra_Vector xs(OperatorRangeMap());
134 Epetra_Vector y(OperatorDomainMap());
135
136 // test operator many times
137 bool result = true;
138 double diffNorm = 0.0, trueNorm = 0.0;
139 for (int i = 0; i < count; i++) {
140 xf.PutScalar(0.0);
141 xs.PutScalar(0.0);
142 y.Random();
143
144 // apply operator
145 Apply(y, xs); // xs = A*y
146 fullContent_->Apply(y, xf); // xf = A*y
147
148 // compute norms
149 xs.Update(-1.0, xf, 1.0);
150 xs.Norm2(&diffNorm);
151 xf.Norm2(&trueNorm);
152
153 // check result
154 result &= (diffNorm / trueNorm < tol);
155 }
156 return result;
157}
158
159} // end namespace Epetra
160} // end namespace Teko
Class that describes how a flat blocked operator should be reordered.
virtual void SetContent(const std::vector< std::vector< int > > &vars, const Teuchos::RCP< const Epetra_Operator > &content)
void Reorder(const BlockReorderManager &brm)
BlockedEpetraOperator(const std::vector< std::vector< int > > &vars, const Teuchos::RCP< const Epetra_Operator > &content, const std::string &label="<ANYM>")
void RemoveReording()
Remove any reordering on this object.
virtual void WriteBlocks(const std::string &prefix) const
bool testAgainstFullOperator(int count, double tol) const
Helps perform sanity checks.
const RCP< const Thyra::LinearOpBase< double > > getThyraOp() const
Return the thyra operator associated with this wrapper.