Teko Version of the Day
Loading...
Searching...
No Matches
Teko_BlockedTpetraOperator.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_BlockedTpetraOperator.hpp"
11#include "Teko_TpetraBlockedMappingStrategy.hpp"
12#include "Teko_TpetraReorderedMappingStrategy.hpp"
13
14#include "Teuchos_VerboseObject.hpp"
15
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"
22
23#include "MatrixMarket_Tpetra.hpp"
24
25#include "Teko_Utilities.hpp"
26
27namespace Teko {
28namespace TpetraHelpers {
29
30using Teuchos::RCP;
31using Teuchos::rcp;
32using Teuchos::rcp_dynamic_cast;
33
35 const std::vector<std::vector<GO> >& vars,
36 const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> >& content, const std::string& label)
37 : Teko::TpetraHelpers::TpetraOperatorWrapper(), label_(label) {
38 SetContent(vars, content);
39}
40
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_);
48
49 // build thyra operator
50 BuildBlockedOperator();
51}
52
53void BlockedTpetraOperator::BuildBlockedOperator() {
54 TEUCHOS_ASSERT(blockedMapping_ != Teuchos::null);
55
56 // get a CRS matrix
57 const RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> > crsContent =
58 rcp_dynamic_cast<const Tpetra::CrsMatrix<ST, LO, GO, NT> >(fullContent_);
59
60 // ask the strategy to build the Thyra operator for you
61 if (blockedOperator_ == Teuchos::null) {
62 blockedOperator_ = blockedMapping_->buildBlockedThyraOp(crsContent, label_);
63 } else {
64 const RCP<Thyra::BlockedLinearOpBase<ST> > blkOp =
65 rcp_dynamic_cast<Thyra::BlockedLinearOpBase<ST> >(blockedOperator_, true);
66 blockedMapping_->rebuildBlockedThyraOp(crsContent, blkOp);
67 }
68
69 // set whatever is returned
70 SetOperator(blockedOperator_, false);
71
72 // reorder if neccessary
73 if (reorderManager_ != Teuchos::null) Reorder(*reorderManager_);
74}
75
76const Teuchos::RCP<const Tpetra::Operator<ST, LO, GO, NT> > BlockedTpetraOperator::GetBlock(
77 int i, int j) const {
78 const RCP<const Thyra::BlockedLinearOpBase<ST> > blkOp =
79 Teuchos::rcp_dynamic_cast<const Thyra::BlockedLinearOpBase<ST> >(getThyraOp());
80
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();
84}
85
90 reorderManager_ = rcp(new BlockReorderManager(brm));
91
92 // build reordered objects
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_);
97
98 RCP<const Thyra::LinearOpBase<ST> > A = buildReorderedLinearOp(*reorderManager_, blockOp);
99
100 // set them as working values
101 SetMapStrategy(reorderMapping);
102 SetOperator(A, false);
103}
104
107 SetMapStrategy(blockedMapping_);
108 SetOperator(blockedOperator_, false);
109 reorderManager_ = Teuchos::null;
110}
111
114void BlockedTpetraOperator::WriteBlocks(const std::string& prefix) const {
115 RCP<Thyra::PhysicallyBlockedLinearOpBase<ST> > blockOp =
116 rcp_dynamic_cast<Thyra::PhysicallyBlockedLinearOpBase<ST> >(blockedOperator_);
117
118 // get size of blocked block operator
119 int rows = Teko::blockRowCount(blockOp);
120
121 for (int i = 0; i < rows; i++) {
122 for (int j = 0; j < rows; j++) {
123 // get the row matrix object
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());
129
130 // write to file
131 Tpetra::MatrixMarket::Writer<Tpetra::CrsMatrix<ST, LO, GO, NT> >::writeSparseFile(
132 formatBlockName(prefix, i, j, rows).c_str(), mat);
133 }
134 }
135}
136
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());
141
142 // test operator many times
143 bool result = true;
144 ST diffNorm = 0.0, trueNorm = 0.0;
145 for (int i = 0; i < count; i++) {
146 xf.putScalar(0.0);
147 xs.putScalar(0.0);
148 y.randomize();
149
150 // apply operator
151 apply(y, xs); // xs = A*y
152 fullContent_->apply(y, xf); // xf = A*y
153
154 // compute norms
155 xs.update(-1.0, xf, 1.0);
156 diffNorm = xs.norm2();
157 trueNorm = xf.norm2();
158
159 // check result
160 result &= (diffNorm / trueNorm < tol);
161 }
162
163 return result;
164}
165
166} // end namespace TpetraHelpers
167} // end namespace Teko
Class that describes how a flat blocked operator should be reordered.
virtual void SetContent(const std::vector< std::vector< GO > > &vars, const Teuchos::RCP< const Tpetra::Operator< ST, LO, GO, NT > > &content)
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.
void RemoveReording()
Remove any reordering on this object.
const RCP< const Thyra::LinearOpBase< ST > > getThyraOp() const
Return the thyra operator associated with this wrapper.