Teko Version of the Day
Loading...
Searching...
No Matches
Teko_InvLSCStrategy.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 "NS/Teko_InvLSCStrategy.hpp"
11
12#include "Thyra_DefaultDiagonalLinearOp.hpp"
13#include "Thyra_VectorStdOps.hpp"
14
15#include "Teko_ConfigDefs.hpp"
16
17#ifdef TEKO_HAVE_EPETRA
18#include "Thyra_EpetraThyraWrappers.hpp"
19#include "Thyra_get_Epetra_Operator.hpp"
20#include "Thyra_EpetraLinearOp.hpp"
21#include "Epetra_Vector.h"
22#include "Epetra_Map.h"
23#include "EpetraExt_RowMatrixOut.h"
24#include "EpetraExt_MultiVectorOut.h"
25#include "EpetraExt_VectorOut.h"
26#include "Teko_EpetraHelpers.hpp"
27#include "Teko_EpetraOperatorWrapper.hpp"
28#endif
29
30#include "Teuchos_Time.hpp"
31#include "Teuchos_TimeMonitor.hpp"
32
33// Teko includes
34#include "Teko_Utilities.hpp"
35#include "NS/Teko_LSCPreconditionerFactory.hpp"
36
37#include "Teko_TpetraHelpers.hpp"
38
39#include "Thyra_TpetraLinearOp.hpp"
40
41using Teuchos::RCP;
42using Teuchos::rcp_const_cast;
43using Teuchos::rcp_dynamic_cast;
44
45namespace Teko {
46namespace NS {
47
49// InvLSCStrategy Implementation
51
52// constructors
54InvLSCStrategy::InvLSCStrategy()
55 : massMatrix_(Teuchos::null),
56 invFactoryF_(Teuchos::null),
57 invFactoryS_(Teuchos::null),
58 eigSolveParam_(5),
59 rowZeroingNeeded_(false),
60 useFullLDU_(false),
61 useMass_(false),
62 useLumping_(false),
63 useWScaling_(false),
64 scaleType_(Diagonal),
65 isSymmetric_(true),
66 assumeStable_(false) {}
67
68InvLSCStrategy::InvLSCStrategy(const Teuchos::RCP<InverseFactory>& factory, bool rzn)
69 : massMatrix_(Teuchos::null),
70 invFactoryF_(factory),
71 invFactoryS_(factory),
72 eigSolveParam_(5),
73 rowZeroingNeeded_(rzn),
74 useFullLDU_(false),
75 useMass_(false),
76 useLumping_(false),
77 useWScaling_(false),
78 scaleType_(Diagonal),
79 isSymmetric_(true),
80 assumeStable_(false) {}
81
82InvLSCStrategy::InvLSCStrategy(const Teuchos::RCP<InverseFactory>& invFactF,
83 const Teuchos::RCP<InverseFactory>& invFactS, bool rzn)
84 : massMatrix_(Teuchos::null),
85 invFactoryF_(invFactF),
86 invFactoryS_(invFactS),
87 eigSolveParam_(5),
88 rowZeroingNeeded_(rzn),
89 useFullLDU_(false),
90 useMass_(false),
91 useLumping_(false),
92 useWScaling_(false),
93 scaleType_(Diagonal),
94 isSymmetric_(true),
95 assumeStable_(false) {}
96
97InvLSCStrategy::InvLSCStrategy(const Teuchos::RCP<InverseFactory>& factory, LinearOp& mass,
98 bool rzn)
99 : massMatrix_(mass),
100 invFactoryF_(factory),
101 invFactoryS_(factory),
102 eigSolveParam_(5),
103 rowZeroingNeeded_(rzn),
104 useFullLDU_(false),
105 useMass_(false),
106 useLumping_(false),
107 useWScaling_(false),
108 scaleType_(Diagonal),
109 isSymmetric_(true),
110 assumeStable_(false) {}
111
112InvLSCStrategy::InvLSCStrategy(const Teuchos::RCP<InverseFactory>& invFactF,
113 const Teuchos::RCP<InverseFactory>& invFactS, LinearOp& mass,
114 bool rzn)
115 : massMatrix_(mass),
116 invFactoryF_(invFactF),
117 invFactoryS_(invFactS),
118 eigSolveParam_(5),
119 rowZeroingNeeded_(rzn),
120 useFullLDU_(false),
121 useMass_(false),
122 useLumping_(false),
123 useWScaling_(false),
124 scaleType_(Diagonal),
125 isSymmetric_(true),
126 assumeStable_(false) {}
127
129
130void InvLSCStrategy::buildState(BlockedLinearOp& A, BlockPreconditionerState& state) const {
131 Teko_DEBUG_SCOPE("InvLSCStrategy::buildState", 10);
132
133 LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
134 TEUCHOS_ASSERT(lscState != 0);
135
136 // if neccessary save state information
137 if (not lscState->isInitialized()) {
138 Teko_DEBUG_EXPR(Teuchos::Time timer(""));
139
140 // construct operators
141 {
142 Teko_DEBUG_SCOPE("LSC::buildState constructing operators", 1);
143 Teko_DEBUG_EXPR(timer.start(true));
144
145 initializeState(A, lscState);
146
147 Teko_DEBUG_EXPR(timer.stop());
148 Teko_DEBUG_MSG("LSC::buildState BuildOpsTime = " << timer.totalElapsedTime(), 1);
149 }
150
151 // Build the inverses
152 {
153 Teko_DEBUG_SCOPE("LSC::buildState calculating inverses", 1);
154 Teko_DEBUG_EXPR(timer.start(true));
155
156 computeInverses(A, lscState);
157
158 Teko_DEBUG_EXPR(timer.stop());
159 Teko_DEBUG_MSG("LSC::buildState BuildInvTime = " << timer.totalElapsedTime(), 1);
160 }
161 }
162}
163
164// functions inherited from LSCStrategy
165LinearOp InvLSCStrategy::getInvBQBt(const BlockedLinearOp& /* A */,
166 BlockPreconditionerState& state) const {
167 return state.getInverse("invBQBtmC");
168}
169
170LinearOp InvLSCStrategy::getInvBHBt(const BlockedLinearOp& /* A */,
171 BlockPreconditionerState& state) const {
172 return state.getInverse("invBHBtmC");
173}
174
175LinearOp InvLSCStrategy::getInvF(const BlockedLinearOp& /* A */,
176 BlockPreconditionerState& state) const {
177 return state.getInverse("invF");
178}
179
180LinearOp InvLSCStrategy::getOuterStabilization(const BlockedLinearOp& /* A */,
181 BlockPreconditionerState& state) const {
182 LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
183 TEUCHOS_ASSERT(lscState != 0);
184 TEUCHOS_ASSERT(lscState->isInitialized())
185
186 return lscState->aiD_;
187}
188
189LinearOp InvLSCStrategy::getInvMass(const BlockedLinearOp& /* A */,
190 BlockPreconditionerState& state) const {
191 LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
192 TEUCHOS_ASSERT(lscState != 0);
193 TEUCHOS_ASSERT(lscState->isInitialized())
194
195 return lscState->invMass_;
196}
197
198LinearOp InvLSCStrategy::getHScaling(const BlockedLinearOp& A,
199 BlockPreconditionerState& state) const {
200 if (hScaling_ != Teuchos::null) return hScaling_;
201 return getInvMass(A, state);
202}
203
205void InvLSCStrategy::initializeState(const BlockedLinearOp& A, LSCPrecondState* state) const {
206 Teko_DEBUG_SCOPE("InvLSCStrategy::initializeState", 10);
207
208 const LinearOp F = getBlock(0, 0, A);
209 const LinearOp Bt = getBlock(0, 1, A);
210 const LinearOp B = getBlock(1, 0, A);
211 const LinearOp C = getBlock(1, 1, A);
212
213 LinearOp D = B;
214 LinearOp G = isSymmetric_ ? Bt : adjoint(D);
215
216 bool isStabilized = assumeStable_ ? false : (not isZeroOp(C));
217
218 // The logic follows like this
219 // if there is no mass matrix available --> build from F
220 // if there is a mass matrix and the inverse hasn't yet been built
221 // --> build from the mass matrix
222 // otherwise, there is already an invMass_ matrix that is appropriate
223 // --> use that one
224 if (massMatrix_ == Teuchos::null) {
225 Teko_DEBUG_MSG(
226 "LSC::initializeState Build Scaling <F> type \"" << getDiagonalName(scaleType_) << "\"", 1);
227 state->invMass_ = getInvDiagonalOp(F, scaleType_);
228 } else if (state->invMass_ == Teuchos::null) {
229 Teko_DEBUG_MSG(
230 "LSC::initializeState Build Scaling <mass> type \"" << getDiagonalName(scaleType_) << "\"",
231 1);
232 state->invMass_ = getInvDiagonalOp(massMatrix_, scaleType_);
233 }
234 // else "invMass_" should be set and there is no reason to rebuild it
235
236 // compute BQBt
237 state->BQBt_ = explicitMultiply(B, state->invMass_, Bt, state->BQBt_);
238 Teko_DEBUG_MSG("Computed BQBt", 10);
239
240 // if there is no H-Scaling
241 if (wScaling_ != Teuchos::null && hScaling_ == Teuchos::null) {
242 // from W vector build H operator scaling
243 RCP<const Thyra::VectorBase<double> > w = wScaling_->col(0);
244 RCP<const Thyra::VectorBase<double> > iQu =
245 rcp_dynamic_cast<const Thyra::DiagonalLinearOpBase<double> >(state->invMass_)->getDiag();
246 RCP<Thyra::VectorBase<double> > h = Thyra::createMember(iQu->space());
247
248 Thyra::put_scalar(0.0, h.ptr());
249 Thyra::ele_wise_prod(1.0, *w, *iQu, h.ptr());
250 hScaling_ = Teuchos::rcp(new Thyra::DefaultDiagonalLinearOp<double>(h));
251 }
252
253 LinearOp H = hScaling_;
254 if (H == Teuchos::null && not isSymmetric_) H = state->invMass_;
255
256 // setup the scaling operator
257 if (H == Teuchos::null)
258 state->BHBt_ = state->BQBt_;
259 else {
260 RCP<Teuchos::Time> time =
261 Teuchos::TimeMonitor::getNewTimer("InvLSCStrategy::initializeState Build BHBt");
262 Teuchos::TimeMonitor timer(*time);
263
264 // compute BHBt
265 state->BHBt_ = explicitMultiply(D, H, G, state->BHBt_);
266 }
267
268 // if this is a stable discretization...we are done!
269 if (not isStabilized) {
270 state->addInverse("BQBtmC", state->BQBt_);
271 state->addInverse("BHBtmC", state->BHBt_);
272 state->gamma_ = 0.0;
273 state->alpha_ = 0.0;
274 state->aiD_ = Teuchos::null;
275
276 state->setInitialized(true);
277
278 return;
279 }
280
281 // for Epetra_CrsMatrix...zero out certain rows: this ensures spectral radius is correct
282 LinearOp modF = F;
283 if (!Teko::TpetraHelpers::isTpetraLinearOp(F)) { // Epetra
284#ifdef TEKO_HAVE_EPETRA
285 const RCP<const Epetra_Operator> epF = Thyra::get_Epetra_Operator(*F);
286 if (epF != Teuchos::null && rowZeroingNeeded_) {
287 // try to get a CRS matrix
288 const RCP<const Epetra_CrsMatrix> crsF = rcp_dynamic_cast<const Epetra_CrsMatrix>(epF);
289
290 // if it is a CRS matrix get rows that need to be zeroed
291 if (crsF != Teuchos::null) {
292 std::vector<int> zeroIndices;
293
294 // get rows in need of zeroing
295 Teko::Epetra::identityRowIndices(crsF->RowMap(), *crsF, zeroIndices);
296
297 // build an operator that zeros those rows
298 modF = Thyra::epetraLinearOp(rcp(new Teko::Epetra::ZeroedOperator(zeroIndices, crsF)));
299 }
300 }
301#else
302 throw std::logic_error(
303 "InvLSCStrategy::initializeState is trying to use "
304 "Epetra code, but TEKO is not built with Epetra!");
305#endif
306 } else { // Tpetra
307 ST scalar = 0.0;
308 bool transp = false;
309 RCP<const Tpetra::CrsMatrix<ST, LO, GO, NT> > crsF =
310 Teko::TpetraHelpers::getTpetraCrsMatrix(F, &scalar, &transp);
311
312 std::vector<GO> zeroIndices;
313
314 // get rows in need of zeroing
315 Teko::TpetraHelpers::identityRowIndices(*crsF->getRowMap(), *crsF, zeroIndices);
316
317 // build an operator that zeros those rows
318 modF = Thyra::tpetraLinearOp<ST, LO, GO, NT>(
319 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(crsF->getDomainMap()),
320 Thyra::tpetraVectorSpace<ST, LO, GO, NT>(crsF->getRangeMap()),
321 rcp(new Teko::TpetraHelpers::ZeroedOperator(zeroIndices, crsF)));
322 }
323
324 // compute gamma
325 Teko_DEBUG_MSG("Calculating gamma", 10);
326 LinearOp iQuF = multiply(state->invMass_, modF);
327
328 // do 6 power iterations to compute spectral radius: EHSST2007 Eq. 4.28
329 Teko::LinearOp stabMatrix; // this is the pressure stabilization matrix to use
330 state->gamma_ = std::fabs(Teko::computeSpectralRad(iQuF, 5e-2, false, eigSolveParam_)) / 3.0;
331 Teko_DEBUG_MSG("Calculated gamma", 10);
332 if (userPresStabMat_ != Teuchos::null) {
333 Teko::LinearOp invDGl = Teko::getInvDiagonalOp(userPresStabMat_);
334 Teko::LinearOp gammaOp = multiply(invDGl, C);
335 state->gamma_ *= std::fabs(Teko::computeSpectralRad(gammaOp, 5e-2, false, eigSolveParam_));
336 stabMatrix = userPresStabMat_;
337 } else
338 stabMatrix = C;
339
340 // compute alpha scaled inv(D): EHSST2007 Eq. 4.29
341 // construct B_idF_Bt and save it for refilling later: This could reuse BQBt graph
342 LinearOp invDiagF = getInvDiagonalOp(F);
343 Teko::ModifiableLinearOp modB_idF_Bt = state->getInverse("BidFBt");
344 modB_idF_Bt = explicitMultiply(B, invDiagF, Bt, modB_idF_Bt);
345 state->addInverse("BidFBt", modB_idF_Bt);
346 const LinearOp B_idF_Bt = modB_idF_Bt;
347
348 MultiVector vec_D = getDiagonal(B_idF_Bt); // this memory could be reused
349 update(-1.0, getDiagonal(C), 1.0, vec_D); // vec_D = diag(B*inv(diag(F))*Bt)-diag(C)
350 const LinearOp invD = buildInvDiagonal(vec_D, "inv(D)");
351
352 Teko_DEBUG_MSG("Calculating alpha", 10);
353 const LinearOp BidFBtidD = multiply<double>(B_idF_Bt, invD);
354 double num = std::fabs(Teko::computeSpectralRad(BidFBtidD, 5e-2, false, eigSolveParam_));
355 Teko_DEBUG_MSG("Calculated alpha", 10);
356 state->alpha_ = 1.0 / num;
357 state->aiD_ = Thyra::scale(state->alpha_, invD);
358
359 // now build B*Q*Bt-gamma*C
360 Teko::ModifiableLinearOp BQBtmC = state->getInverse("BQBtmC");
361 BQBtmC = explicitAdd(state->BQBt_, scale(-state->gamma_, stabMatrix), BQBtmC);
362 state->addInverse("BQBtmC", BQBtmC);
363
364 // now build B*H*Bt-gamma*C
365 Teko::ModifiableLinearOp BHBtmC = state->getInverse("BHBtmC");
366 if (H == Teuchos::null)
367 BHBtmC = BQBtmC;
368 else {
369 BHBtmC = explicitAdd(state->BHBt_, scale(-state->gamma_, stabMatrix), BHBtmC);
370 }
371 state->addInverse("BHBtmC", BHBtmC);
372
373 Teko_DEBUG_MSG_BEGIN(5) DEBUG_STREAM << "LSC Gamma Parameter = " << state->gamma_ << std::endl;
374 DEBUG_STREAM << "LSC Alpha Parameter = " << state->alpha_ << std::endl;
375 Teko_DEBUG_MSG_END()
376
377 state->setInitialized(true);
378}
379
385void InvLSCStrategy::computeInverses(const BlockedLinearOp& A, LSCPrecondState* state) const {
386 Teko_DEBUG_SCOPE("InvLSCStrategy::computeInverses", 10);
387 Teko_DEBUG_EXPR(Teuchos::Time invTimer(""));
388
389 const LinearOp F = getBlock(0, 0, A);
390
392
393 // (re)build the inverse of F
394 Teko_DEBUG_MSG("LSC::computeInverses Building inv(F)", 1);
395 Teko_DEBUG_EXPR(invTimer.start(true));
396 InverseLinearOp invF = state->getInverse("invF");
397 if (invF == Teuchos::null) {
398 invF = buildInverse(*invFactoryF_, F);
399 state->addInverse("invF", invF);
400 } else {
401 rebuildInverse(*invFactoryF_, F, invF);
402 }
403 Teko_DEBUG_EXPR(invTimer.stop());
404 Teko_DEBUG_MSG("LSC::computeInverses GetInvF = " << invTimer.totalElapsedTime(), 1);
405
407
408 // (re)build the inverse of BQBt
409 Teko_DEBUG_MSG("LSC::computeInverses Building inv(BQBtmC)", 1);
410 Teko_DEBUG_EXPR(invTimer.start(true));
411 const LinearOp BQBt = state->getInverse("BQBtmC");
412 InverseLinearOp invBQBt = state->getInverse("invBQBtmC");
413 if (invBQBt == Teuchos::null) {
414 invBQBt = buildInverse(*invFactoryS_, BQBt);
415 state->addInverse("invBQBtmC", invBQBt);
416 } else {
417 rebuildInverse(*invFactoryS_, BQBt, invBQBt);
418 }
419 Teko_DEBUG_EXPR(invTimer.stop());
420 Teko_DEBUG_MSG("LSC::computeInverses GetInvBQBt = " << invTimer.totalElapsedTime(), 1);
421
423
424 // Compute the inverse of BHBt or just use BQBt
425 ModifiableLinearOp invBHBt = state->getInverse("invBHBtmC");
426 if (hScaling_ != Teuchos::null || not isSymmetric_) {
427 // (re)build the inverse of BHBt
428 Teko_DEBUG_MSG("LSC::computeInverses Building inv(BHBtmC)", 1);
429 Teko_DEBUG_EXPR(invTimer.start(true));
430 const LinearOp BHBt = state->getInverse("BHBtmC");
431 if (invBHBt == Teuchos::null) {
432 invBHBt = buildInverse(*invFactoryS_, BHBt);
433 state->addInverse("invBHBtmC", invBHBt);
434 } else {
435 rebuildInverse(*invFactoryS_, BHBt, invBHBt);
436 }
437 Teko_DEBUG_EXPR(invTimer.stop());
438 Teko_DEBUG_MSG("LSC::computeInverses GetInvBHBt = " << invTimer.totalElapsedTime(), 1);
439 } else if (invBHBt == Teuchos::null) {
440 // just use the Q version
441 state->addInverse("invBHBtmC", invBQBt);
442 }
443}
444
446void InvLSCStrategy::initializeFromParameterList(const Teuchos::ParameterList& pl,
447 const InverseLibrary& invLib) {
448 // get string specifying inverse
449 std::string invStr = "", invVStr = "", invPStr = "";
450 bool rowZeroing = true;
451 bool useLDU = false;
452 scaleType_ = Diagonal;
453
454 // "parse" the parameter list
455 if (pl.isParameter("Inverse Type")) invStr = pl.get<std::string>("Inverse Type");
456 if (pl.isParameter("Inverse Velocity Type"))
457 invVStr = pl.get<std::string>("Inverse Velocity Type");
458 if (pl.isParameter("Inverse Pressure Type"))
459 invPStr = pl.get<std::string>("Inverse Pressure Type");
460 if (pl.isParameter("Ignore Boundary Rows")) rowZeroing = pl.get<bool>("Ignore Boundary Rows");
461 if (pl.isParameter("Use LDU")) useLDU = pl.get<bool>("Use LDU");
462 if (pl.isParameter("Use Mass Scaling")) useMass_ = pl.get<bool>("Use Mass Scaling");
463 // if(pl.isParameter("Use Lumping"))
464 // useLumping_ = pl.get<bool>("Use Lumping");
465 if (pl.isParameter("Use W-Scaling")) useWScaling_ = pl.get<bool>("Use W-Scaling");
466 if (pl.isParameter("Eigen Solver Iterations"))
467 eigSolveParam_ = pl.get<int>("Eigen Solver Iterations");
468 if (pl.isParameter("Scaling Type")) {
469 scaleType_ = getDiagonalType(pl.get<std::string>("Scaling Type"));
470 TEUCHOS_TEST_FOR_EXCEPT(scaleType_ == NotDiag);
471 }
472 if (pl.isParameter("Assume Stable Discretization"))
473 assumeStable_ = pl.get<bool>("Assume Stable Discretization");
474
475 Teko_DEBUG_MSG_BEGIN(5) DEBUG_STREAM << "LSC Inverse Strategy Parameters: " << std::endl;
476 DEBUG_STREAM << " inv type = \"" << invStr << "\"" << std::endl;
477 DEBUG_STREAM << " inv v type = \"" << invVStr << "\"" << std::endl;
478 DEBUG_STREAM << " inv p type = \"" << invPStr << "\"" << std::endl;
479 DEBUG_STREAM << " bndry rows = " << rowZeroing << std::endl;
480 DEBUG_STREAM << " use ldu = " << useLDU << std::endl;
481 DEBUG_STREAM << " use mass = " << useMass_ << std::endl;
482 DEBUG_STREAM << " use w-scaling = " << useWScaling_ << std::endl;
483 DEBUG_STREAM << " assume stable = " << assumeStable_ << std::endl;
484 DEBUG_STREAM << " scale type = " << getDiagonalName(scaleType_) << std::endl;
485 DEBUG_STREAM << "LSC Inverse Strategy Parameter list: " << std::endl;
486 pl.print(DEBUG_STREAM);
487 Teko_DEBUG_MSG_END()
488
489 // set defaults as needed
490#if defined(Teko_ENABLE_Amesos)
491 if (invStr == "") invStr = "Amesos";
492#elif defined(Teko_ENABLE_Amesos2)
493 if (invStr == "") invStr = "Amesos2";
494#endif
495 if (invVStr == "") invVStr = invStr;
496 if (invPStr == "") invPStr = invStr;
497
498 // build velocity inverse factory
499 invFactoryF_ = invLib.getInverseFactory(invVStr);
500 invFactoryS_ = invFactoryF_; // by default these are the same
501 if (invVStr != invPStr) // if different, build pressure inverse factory
502 invFactoryS_ = invLib.getInverseFactory(invPStr);
503
504 // set other parameters
505 setUseFullLDU(useLDU);
506 setRowZeroing(rowZeroing);
507
508 if (useMass_) {
509 Teuchos::RCP<Teko::RequestHandler> rh = getRequestHandler();
510 rh->preRequest<Teko::LinearOp>(Teko::RequestMesg("Velocity Mass Matrix"));
511 Teko::LinearOp mass = rh->request<Teko::LinearOp>(Teko::RequestMesg("Velocity Mass Matrix"));
512 setMassMatrix(mass);
513 }
514}
515
517Teuchos::RCP<Teuchos::ParameterList> InvLSCStrategy::getRequestedParameters() const {
518 Teuchos::RCP<Teuchos::ParameterList> result;
519 Teuchos::RCP<Teuchos::ParameterList> pl = rcp(new Teuchos::ParameterList());
520
521 // grab parameters from F solver
522 RCP<Teuchos::ParameterList> fList = invFactoryF_->getRequestedParameters();
523 if (fList != Teuchos::null) {
524 Teuchos::ParameterList::ConstIterator itr;
525 for (itr = fList->begin(); itr != fList->end(); ++itr) pl->setEntry(itr->first, itr->second);
526 result = pl;
527 }
528
529 // grab parameters from S solver
530 RCP<Teuchos::ParameterList> sList = invFactoryS_->getRequestedParameters();
531 if (sList != Teuchos::null) {
532 Teuchos::ParameterList::ConstIterator itr;
533 for (itr = sList->begin(); itr != sList->end(); ++itr) pl->setEntry(itr->first, itr->second);
534 result = pl;
535 }
536
537 // use the mass matrix
538 if (useWScaling_) {
539 pl->set<Teko::LinearOp>("W-Scaling Vector", Teuchos::null, "W-Scaling Vector");
540 result = pl;
541 }
542
543 return result;
544}
545
547bool InvLSCStrategy::updateRequestedParameters(const Teuchos::ParameterList& pl) {
548 Teko_DEBUG_SCOPE("InvLSCStrategy::updateRequestedParameters", 10);
549 bool result = true;
550
551 // update requested parameters in solvers
552 result &= invFactoryF_->updateRequestedParameters(pl);
553 result &= invFactoryS_->updateRequestedParameters(pl);
554
555 // use W scaling matrix
556 if (useWScaling_) {
557 Teko::MultiVector wScale = pl.get<Teko::MultiVector>("W-Scaling Vector");
558
559 if (wScale == Teuchos::null)
560 result &= false;
561 else
562 setWScaling(wScale);
563 }
564
565 return result;
566}
567
568} // end namespace NS
569} // 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.
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl, const InverseLibrary &invLib)
Initialize from a parameter list.
virtual LinearOp getInvMass(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setMassMatrix(const LinearOp &mass)
set the mass matrix to use in computing the scaling
virtual LinearOp getHScaling(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setUseFullLDU(bool val)
Set to true to use the Full LDU decomposition, false otherwise.
virtual LinearOp getInvF(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setWScaling(const MultiVector &wScaling)
virtual void buildState(BlockedLinearOp &A, BlockPreconditionerState &state) const
Functions inherited from LSCStrategy.
virtual LinearOp getOuterStabilization(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual bool updateRequestedParameters(const Teuchos::ParameterList &pl)
For assiting in construction of the preconditioner.
virtual void initializeState(const BlockedLinearOp &A, LSCPrecondState *state) const
Initialize the state object using this blocked linear operator.
virtual LinearOp getInvBHBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual Teuchos::RCP< Teuchos::ParameterList > getRequestedParameters() const
For assiting in construction of the preconditioner.
void computeInverses(const BlockedLinearOp &A, LSCPrecondState *state) const
virtual void setRowZeroing(bool val)
Set to true to zero the rows of F when computing the spectral radius.
virtual LinearOp getInvBQBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
Preconditioner state for the LSC factory.
LinearOp invMass_
Inverse mass operator ( ).
Teuchos::RCP< RequestHandler > getRequestHandler() const
This method gets the request handler uses by this object.
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.