Teko Version of the Day
Loading...
Searching...
No Matches
Teko_StridedEpetraOperator.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 "Epetra/Teko_StridedEpetraOperator.hpp"
11#include "Epetra/Teko_StridedMappingStrategy.hpp"
12#include "Epetra/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 "Epetra_Vector.h"
25#include "EpetraExt_MultiVectorOut.h"
26#include "EpetraExt_RowMatrixOut.h"
27
28#include "Teko_Utilities.hpp"
29
30namespace Teko {
31namespace Epetra {
32
33using Teuchos::RCP;
34using Teuchos::rcp;
35using Teuchos::rcp_dynamic_cast;
36
37StridedEpetraOperator::StridedEpetraOperator(int numVars,
38 const Teuchos::RCP<const Epetra_Operator>& content,
39 const std::string& label)
40 : Teko::Epetra::EpetraOperatorWrapper(), label_(label) {
41 std::vector<int> vars;
42
43 // build vector describing the sub maps
44 for (int i = 0; i < numVars; i++) vars.push_back(1);
45
46 SetContent(vars, content);
47}
48
49StridedEpetraOperator::StridedEpetraOperator(const std::vector<int>& vars,
50 const Teuchos::RCP<const Epetra_Operator>& content,
51 const std::string& label)
52 : Teko::Epetra::EpetraOperatorWrapper(), label_(label) {
53 SetContent(vars, content);
54}
55
56void StridedEpetraOperator::SetContent(const std::vector<int>& vars,
57 const Teuchos::RCP<const Epetra_Operator>& content) {
58 fullContent_ = content;
59 stridedMapping_ = rcp(new StridedMappingStrategy(
60 vars, Teuchos::rcpFromRef(fullContent_->OperatorDomainMap()), fullContent_->Comm()));
61 SetMapStrategy(stridedMapping_);
62
63 // build thyra operator
64 BuildBlockedOperator();
65}
66
67void StridedEpetraOperator::BuildBlockedOperator() {
68 TEUCHOS_ASSERT(stridedMapping_ != Teuchos::null);
69
70 // get a CRS matrix
71 const RCP<const Epetra_CrsMatrix> crsContent =
72 rcp_dynamic_cast<const Epetra_CrsMatrix>(fullContent_);
73
74 // ask the strategy to build the Thyra operator for you
75 if (stridedOperator_ == Teuchos::null) {
76 stridedOperator_ = stridedMapping_->buildBlockedThyraOp(crsContent, label_);
77 } else {
78 const RCP<Thyra::BlockedLinearOpBase<double> > blkOp =
79 rcp_dynamic_cast<Thyra::BlockedLinearOpBase<double> >(stridedOperator_, true);
80 stridedMapping_->rebuildBlockedThyraOp(crsContent, blkOp);
81 }
82
83 // set whatever is returned
84 SetOperator(stridedOperator_, false);
85
86 // reorder if neccessary
87 if (reorderManager_ != Teuchos::null) Reorder(*reorderManager_);
88}
89
90const Teuchos::RCP<const Epetra_Operator> StridedEpetraOperator::GetBlock(int i, int j) const {
91 const RCP<const Thyra::BlockedLinearOpBase<double> > blkOp =
92 Teuchos::rcp_dynamic_cast<const Thyra::BlockedLinearOpBase<double> >(getThyraOp());
93
94 return Thyra::get_Epetra_Operator(*blkOp->getBlock(i, j));
95}
96
100void StridedEpetraOperator::Reorder(const BlockReorderManager& brm) {
101 reorderManager_ = rcp(new BlockReorderManager(brm));
102
103 // build reordered objects
104 RCP<const MappingStrategy> reorderMapping =
105 rcp(new ReorderedMappingStrategy(*reorderManager_, stridedMapping_));
106 RCP<const Thyra::BlockedLinearOpBase<double> > blockOp =
107 rcp_dynamic_cast<const Thyra::BlockedLinearOpBase<double> >(stridedOperator_);
108
109 RCP<const Thyra::LinearOpBase<double> > A = buildReorderedLinearOp(*reorderManager_, blockOp);
110
111 // set them as working values
112 SetMapStrategy(reorderMapping);
113 SetOperator(A, false);
114}
115
117void StridedEpetraOperator::RemoveReording() {
118 SetMapStrategy(stridedMapping_);
119 SetOperator(stridedOperator_, false);
120 reorderManager_ = Teuchos::null;
121}
122
125void StridedEpetraOperator::WriteBlocks(const std::string& prefix) const {
126 RCP<Thyra::PhysicallyBlockedLinearOpBase<double> > blockOp =
127 rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<double> >(stridedOperator_);
128
129 // get size of strided block operator
130 int rows = Teko::blockRowCount(blockOp);
131
132 for (int i = 0; i < rows; i++) {
133 for (int j = 0; j < rows; j++) {
134 // get the row matrix object (Note: can't use "GetBlock" method b/c matrix might be reordered)
135 RCP<const Epetra_RowMatrix> mat = Teuchos::rcp_dynamic_cast<const Epetra_RowMatrix>(
136 Thyra::get_Epetra_Operator(*blockOp->getBlock(i, j)));
137
138 // write to file
139 EpetraExt::RowMatrixToMatrixMarketFile(formatBlockName(prefix, i, j, rows).c_str(), *mat);
140 }
141 }
142}
143
151std::string StridedEpetraOperator::PrintNorm(const eNormType& nrmType, const char newline) {
152 BlockedLinearOp blockOp = toBlockedLinearOp(stridedOperator_);
153
154 // get size of strided block operator
155 int rowCount = Teko::blockRowCount(blockOp);
156 int colCount = Teko::blockRowCount(blockOp);
157
158 std::stringstream ss;
159 ss.precision(4);
160 ss << std::scientific;
161 for (int row = 0; row < rowCount; row++) {
162 for (int col = 0; col < colCount; col++) {
163 // get the row matrix object (Note: can't use "GetBlock" method b/c matrix might be reordered)
164 RCP<const Epetra_CrsMatrix> mat = Teuchos::rcp_dynamic_cast<const Epetra_CrsMatrix>(
165 Thyra::get_Epetra_Operator(*blockOp->getBlock(row, col)));
166
167 // compute the norm
168 double norm = 0.0;
169 switch (nrmType) {
170 case Inf: norm = mat->NormInf(); break;
171 case One: norm = mat->NormOne(); break;
172 case Frobenius: norm = mat->NormFrobenius(); break;
173 default:
174 TEUCHOS_TEST_FOR_EXCEPTION(true, std::runtime_error,
175 "StridedEpetraOperator::eNormType incorrectly specified");
176 }
177
178 ss << norm << " ";
179 }
180 ss << newline;
181 }
182
183 return ss.str();
184}
185
186bool StridedEpetraOperator::testAgainstFullOperator(int count, double tol) const {
187 Epetra_Vector xf(OperatorRangeMap());
188 Epetra_Vector xs(OperatorRangeMap());
189 Epetra_Vector y(OperatorDomainMap());
190
191 // test operator many times
192 bool result = true;
193 double diffNorm = 0.0, trueNorm = 0.0;
194 for (int i = 0; i < count; i++) {
195 xf.PutScalar(0.0);
196 xs.PutScalar(0.0);
197 y.Random();
198
199 // apply operator
200 Apply(y, xs); // xs = A*y
201 fullContent_->Apply(y, xf); // xf = A*y
202
203 // compute norms
204 xs.Update(-1.0, xf, 1.0);
205 xs.Norm2(&diffNorm);
206 xf.Norm2(&trueNorm);
207
208 // check result
209 result &= (diffNorm / trueNorm < tol);
210 }
211
212 return result;
213}
214
215} // end namespace Epetra
216} // end namespace Teko
Implements the Epetra_Operator interface with a Thyra LinearOperator. This enables the use of absrtac...