Teko Version of the Day
Loading...
Searching...
No Matches
Teko_LSCSIMPLECStrategy.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_LSCSIMPLECStrategy.hpp"
11
12#include "Thyra_DefaultDiagonalLinearOp.hpp"
13#include "Thyra_VectorStdOps.hpp"
14
15#include "Teuchos_Time.hpp"
16#include "Teuchos_TimeMonitor.hpp"
17
18// Teko includes
19#include "Teko_Utilities.hpp"
20#include "Teko_LSCPreconditionerFactory.hpp"
21
22using Teuchos::RCP;
23using Teuchos::rcp_const_cast;
24using Teuchos::rcp_dynamic_cast;
25
26namespace Teko {
27namespace NS {
28
30// LSCSIMPLECStrategy Implementation
32
33// constructors
35LSCSIMPLECStrategy::LSCSIMPLECStrategy()
36 : invFactoryF_(Teuchos::null),
37 invFactoryS_(Teuchos::null),
38 useFullLDU_(false),
39 scaleType_(Diagonal) {}
40
42
43void LSCSIMPLECStrategy::buildState(BlockedLinearOp& A, BlockPreconditionerState& state) const {
44 Teko_DEBUG_SCOPE("LSCSIMPLECStrategy::buildState", 10);
45
46 LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
47 TEUCHOS_ASSERT(lscState != 0);
48
49 // if neccessary save state information
50 if (not lscState->isInitialized()) {
51 Teko_DEBUG_EXPR(Teuchos::Time timer(""));
52
53 // construct operators
54 {
55 Teko_DEBUG_SCOPE("LSC-SIMPLEC::buildState constructing operators", 1);
56 Teko_DEBUG_EXPR(timer.start(true));
57
58 initializeState(A, lscState);
59
60 Teko_DEBUG_EXPR(timer.stop());
61 Teko_DEBUG_MSG("LSC-SIMPLEC::buildState BuildOpsTime = " << timer.totalElapsedTime(), 1);
62 }
63
64 // Build the inverses
65 {
66 Teko_DEBUG_SCOPE("LSC-SIMPLEC::buildState calculating inverses", 1);
67 Teko_DEBUG_EXPR(timer.start(true));
68
69 computeInverses(A, lscState);
70
71 Teko_DEBUG_EXPR(timer.stop());
72 Teko_DEBUG_MSG("LSC-SIMPLEC::buildState BuildInvTime = " << timer.totalElapsedTime(), 1);
73 }
74 }
75}
76
77// functions inherited from LSCStrategy
78LinearOp LSCSIMPLECStrategy::getInvBQBt(const BlockedLinearOp& /* A */,
79 BlockPreconditionerState& state) const {
80 return state.getInverse("invBQBtmC");
81}
82
83LinearOp LSCSIMPLECStrategy::getInvBHBt(const BlockedLinearOp& /* A */,
84 BlockPreconditionerState& state) const {
85 return state.getInverse("invBQBtmC").getConst();
86}
87
88LinearOp LSCSIMPLECStrategy::getInvF(const BlockedLinearOp& /* A */,
89 BlockPreconditionerState& state) const {
90 return state.getInverse("invF");
91}
92
93LinearOp LSCSIMPLECStrategy::getInnerStabilization(const BlockedLinearOp& A,
94 BlockPreconditionerState& state) const {
95 LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
96 TEUCHOS_ASSERT(lscState != 0);
97 TEUCHOS_ASSERT(lscState->isInitialized())
98
99 const LinearOp C = getBlock(1, 1, A);
100 return scale(-1.0, C);
101}
102
103LinearOp LSCSIMPLECStrategy::getInvMass(const BlockedLinearOp& /* A */,
104 BlockPreconditionerState& state) const {
105 LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
106 TEUCHOS_ASSERT(lscState != 0);
107 TEUCHOS_ASSERT(lscState->isInitialized())
108
109 return lscState->invMass_;
110}
111
112LinearOp LSCSIMPLECStrategy::getHScaling(const BlockedLinearOp& A,
113 BlockPreconditionerState& state) const {
114 return getInvMass(A, state);
115}
116
118void LSCSIMPLECStrategy::initializeState(const BlockedLinearOp& A, LSCPrecondState* state) const {
119 Teko_DEBUG_SCOPE("LSCSIMPLECStrategy::initializeState", 10);
120
121 const LinearOp F = getBlock(0, 0, A);
122 const LinearOp Bt = getBlock(0, 1, A);
123 const LinearOp B = getBlock(1, 0, A);
124 const LinearOp C = getBlock(1, 1, A);
125
126 bool isStabilized = (not isZeroOp(C));
127
128 state->invMass_ = getInvDiagonalOp(F, scaleType_);
129
130 // compute BQBt
131 state->BQBt_ = explicitMultiply(B, state->invMass_, Bt, state->BQBt_);
132 if (isStabilized) {
133 // now build B*Q*Bt-C
134 Teko::ModifiableLinearOp BQBtmC = state->getInverse("BQBtmC");
135 BQBtmC = explicitAdd(state->BQBt_, scale(-1.0, C), BQBtmC);
136 state->addInverse("BQBtmC", BQBtmC);
137 }
138 Teko_DEBUG_MSG("Computed BQBt", 10);
139
140 state->setInitialized(true);
141}
142
148void LSCSIMPLECStrategy::computeInverses(const BlockedLinearOp& A, LSCPrecondState* state) const {
149 Teko_DEBUG_SCOPE("LSCSIMPLECStrategy::computeInverses", 10);
150 Teko_DEBUG_EXPR(Teuchos::Time invTimer(""));
151
152 const LinearOp F = getBlock(0, 0, A);
153
155
156 // (re)build the inverse of F
157 Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses Building inv(F)", 1);
158 Teko_DEBUG_EXPR(invTimer.start(true));
159 InverseLinearOp invF = state->getInverse("invF");
160 if (invF == Teuchos::null) {
161 invF = buildInverse(*invFactoryF_, F);
162 state->addInverse("invF", invF);
163 } else {
164 rebuildInverse(*invFactoryF_, F, invF);
165 }
166 Teko_DEBUG_EXPR(invTimer.stop());
167 Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses GetInvF = " << invTimer.totalElapsedTime(), 1);
168
170
171 // (re)build the inverse of BQBt
172 Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses Building inv(BQBtmC)", 1);
173 Teko_DEBUG_EXPR(invTimer.start(true));
174 const LinearOp BQBt = state->getInverse("BQBtmC");
175 InverseLinearOp invBQBt = state->getInverse("invBQBtmC");
176 if (invBQBt == Teuchos::null) {
177 invBQBt = buildInverse(*invFactoryS_, BQBt);
178 state->addInverse("invBQBtmC", invBQBt);
179 } else {
180 rebuildInverse(*invFactoryS_, BQBt, invBQBt);
181 }
182 Teko_DEBUG_EXPR(invTimer.stop());
183 Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses GetInvBQBt = " << invTimer.totalElapsedTime(), 1);
184}
185
187void LSCSIMPLECStrategy::initializeFromParameterList(const Teuchos::ParameterList& pl,
188 const InverseLibrary& invLib) {
189 // get string specifying inverse
190 std::string invStr = "", invVStr = "", invPStr = "";
191 bool useLDU = false;
192 scaleType_ = Diagonal;
193
194 // "parse" the parameter list
195 if (pl.isParameter("Inverse Type")) invStr = pl.get<std::string>("Inverse Type");
196 if (pl.isParameter("Inverse Velocity Type"))
197 invVStr = pl.get<std::string>("Inverse Velocity Type");
198 if (pl.isParameter("Inverse Pressure Type"))
199 invPStr = pl.get<std::string>("Inverse Pressure Type");
200 if (pl.isParameter("Use LDU")) useLDU = pl.get<bool>("Use LDU");
201 if (pl.isParameter("Scaling Type")) {
202 scaleType_ = getDiagonalType(pl.get<std::string>("Scaling Type"));
203 TEUCHOS_TEST_FOR_EXCEPT(scaleType_ == NotDiag);
204 }
205
206 Teko_DEBUG_MSG_BEGIN(0) DEBUG_STREAM << "LSC Inverse Strategy Parameters: " << std::endl;
207 DEBUG_STREAM << " inv type = \"" << invStr << "\"" << std::endl;
208 DEBUG_STREAM << " inv v type = \"" << invVStr << "\"" << std::endl;
209 DEBUG_STREAM << " inv p type = \"" << invPStr << "\"" << std::endl;
210 DEBUG_STREAM << " use ldu = " << useLDU << std::endl;
211 DEBUG_STREAM << " scale type = " << getDiagonalName(scaleType_) << std::endl;
212 DEBUG_STREAM << "LSC Inverse Strategy Parameter list: " << std::endl;
213 pl.print(DEBUG_STREAM);
214 Teko_DEBUG_MSG_END()
215
216 // set defaults as needed
217#if defined(Teko_ENABLE_Amesos)
218 if (invStr == "") invStr = "Amesos";
219#elif defined(Teko_ENABLE_Amesos2)
220 if (invStr == "") invStr = "Amesos2";
221#endif
222 if (invVStr == "") invVStr = invStr;
223 if (invPStr == "") invPStr = invStr;
224
225 // build velocity inverse factory
226 invFactoryF_ = invLib.getInverseFactory(invVStr);
227 invFactoryS_ = invFactoryF_; // by default these are the same
228 if (invVStr != invPStr) // if different, build pressure inverse factory
229 invFactoryS_ = invLib.getInverseFactory(invPStr);
230
231 // set other parameters
232 setUseFullLDU(useLDU);
233}
234
235} // end namespace NS
236} // end namespace Teko
@ NotDiag
For user convenience, if Teko recieves this value, exceptions will be thrown.
@ Diagonal
Specifies that just the diagonal is used.
MultiVector getBlock(int i, const BlockedMultiVector &bmv)
Get the ith block from a BlockedMultiVector object.
An implementation of a state object for block preconditioners.
Preconditioner state for the LSC factory.
LinearOp invMass_
Inverse mass operator ( ).
virtual void buildState(BlockedLinearOp &A, BlockPreconditionerState &state) const
Functions inherited from LSCStrategy.
virtual LinearOp getInvBHBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual LinearOp getInvF(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl, const InverseLibrary &invLib)
Initialize from a parameter list.
virtual LinearOp getHScaling(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual LinearOp getInvBQBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setUseFullLDU(bool val)
Set to true to use the Full LDU decomposition, false otherwise.
void computeInverses(const BlockedLinearOp &A, LSCPrecondState *state) const
virtual void initializeState(const BlockedLinearOp &A, LSCPrecondState *state) const
Initialize the state object using this blocked linear operator.
virtual LinearOp getInvMass(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void addInverse(const std::string &name, const Teko::InverseLinearOp &ilo)
Add a named inverse to the state object.
virtual void setInitialized(bool init=true)
virtual Teko::InverseLinearOp getInverse(const std::string &name) const
Get a named inverse from the state object.