Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_PresLaplaceLSCStrategy.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_PresLaplaceLSCStrategy.hpp"
11 
12 #include "Thyra_DefaultDiagonalLinearOp.hpp"
13 
14 #include "Teuchos_Time.hpp"
15 #include "Teuchos_TimeMonitor.hpp"
16 
17 // Teko includes
18 #include "Teko_Utilities.hpp"
19 #include "NS/Teko_LSCPreconditionerFactory.hpp"
20 
21 using Teuchos::RCP;
22 using Teuchos::rcp_const_cast;
23 using Teuchos::rcp_dynamic_cast;
24 
25 namespace Teko {
26 namespace NS {
27 
29 // PresLaplaceLSCStrategy Implementation
31 
32 // constructors
34 PresLaplaceLSCStrategy::PresLaplaceLSCStrategy()
35  : invFactoryV_(Teuchos::null),
36  invFactoryP_(Teuchos::null),
37  eigSolveParam_(5),
38  useFullLDU_(false),
39  useMass_(false),
40  scaleType_(AbsRowSum) {}
41 
42 PresLaplaceLSCStrategy::PresLaplaceLSCStrategy(const Teuchos::RCP<InverseFactory>& factory)
43  : invFactoryV_(factory),
44  invFactoryP_(factory),
45  eigSolveParam_(5),
46  useFullLDU_(false),
47  useMass_(false),
48  scaleType_(AbsRowSum) {}
49 
50 PresLaplaceLSCStrategy::PresLaplaceLSCStrategy(const Teuchos::RCP<InverseFactory>& invFactF,
51  const Teuchos::RCP<InverseFactory>& invFactS)
52  : invFactoryV_(invFactF),
53  invFactoryP_(invFactS),
54  eigSolveParam_(5),
55  useFullLDU_(false),
56  useMass_(false),
57  scaleType_(AbsRowSum) {}
58 
60 
61 void PresLaplaceLSCStrategy::buildState(BlockedLinearOp& A, BlockPreconditionerState& state) const {
62  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::buildState", 10);
63 
64  LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
65  TEUCHOS_ASSERT(lscState != 0);
66 
67  // if neccessary save state information
68  if (not lscState->isInitialized()) {
69  Teko_DEBUG_EXPR(Teuchos::Time timer(""));
70 
71  // construct operators
72  {
73  Teko_DEBUG_SCOPE("PL-LSC::buildState constructing operators", 1);
74  Teko_DEBUG_EXPR(timer.start(true));
75 
76  initializeState(A, lscState);
77 
78  Teko_DEBUG_EXPR(timer.stop());
79  Teko_DEBUG_MSG("PL-LSC::buildState BuildOpsTime = " << timer.totalElapsedTime(), 1);
80  }
81 
82  // Build the inverses
83  {
84  Teko_DEBUG_SCOPE("PL-LSC::buildState calculating inverses", 1);
85  Teko_DEBUG_EXPR(timer.start(true));
86 
87  computeInverses(A, lscState);
88 
89  Teko_DEBUG_EXPR(timer.stop());
90  Teko_DEBUG_MSG("PL-LSC::buildState BuildInvTime = " << timer.totalElapsedTime(), 1);
91  }
92  }
93 }
94 
95 // functions inherited from LSCStrategy
96 LinearOp PresLaplaceLSCStrategy::getInvBQBt(const BlockedLinearOp& /* A */,
97  BlockPreconditionerState& state) const {
98  return state.getModifiableOp("invPresLap");
99 }
100 
101 LinearOp PresLaplaceLSCStrategy::getInvBHBt(const BlockedLinearOp& /* A */,
102  BlockPreconditionerState& state) const {
103  return state.getModifiableOp("invPresLap");
104 }
105 
106 LinearOp PresLaplaceLSCStrategy::getInvF(const BlockedLinearOp& /* A */,
107  BlockPreconditionerState& state) const {
108  return state.getModifiableOp("invF");
109 }
110 
111 // LinearOp PresLaplaceLSCStrategy::getInvAlphaD(const BlockedLinearOp & A,BlockPreconditionerState
112 // & state) const
113 LinearOp PresLaplaceLSCStrategy::getOuterStabilization(const BlockedLinearOp& /* A */,
114  BlockPreconditionerState& state) const {
115  LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
116  TEUCHOS_ASSERT(lscState != 0);
117  TEUCHOS_ASSERT(lscState->isInitialized())
118 
119  return lscState->aiD_;
120 }
121 
122 LinearOp PresLaplaceLSCStrategy::getInvMass(const BlockedLinearOp& /* A */,
123  BlockPreconditionerState& state) const {
124  LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
125  TEUCHOS_ASSERT(lscState != 0);
126  TEUCHOS_ASSERT(lscState->isInitialized())
127 
128  return lscState->invMass_;
129 }
130 
131 LinearOp PresLaplaceLSCStrategy::getHScaling(const BlockedLinearOp& A,
132  BlockPreconditionerState& state) const {
133  return getInvMass(A, state);
134 }
135 
137 void PresLaplaceLSCStrategy::initializeState(const BlockedLinearOp& A,
138  LSCPrecondState* state) const {
139  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::initializeState", 10);
140 
141  std::string velMassStr = getVelocityMassString();
142 
143  const LinearOp F = getBlock(0, 0, A);
144  const LinearOp Bt = getBlock(0, 1, A);
145  const LinearOp B = getBlock(1, 0, A);
146  const LinearOp C = getBlock(1, 1, A);
147 
148  LinearOp D = B;
149  LinearOp G = Bt;
150 
151  bool isStabilized = (not isZeroOp(C));
152 
153  // grab operators from state object
154  LinearOp massMatrix = state->getLinearOp(velMassStr);
155 
156  // The logic follows like this
157  // if there is no mass matrix available --> build from F
158  // if there is a mass matrix and the inverse hasn't yet been built
159  // --> build from the mass matrix
160  // otherwise, there is already an invMass_ matrix that is appropriate
161  // --> use that one
162  if (massMatrix == Teuchos::null) {
163  Teko_DEBUG_MSG(
164  "PL-LSC::initializeState Build Scaling <F> type \"" << getDiagonalName(scaleType_) << "\"",
165  1);
166  state->invMass_ = getInvDiagonalOp(F, scaleType_);
167  } else if (state->invMass_ == Teuchos::null) {
168  Teko_DEBUG_MSG("PL-LSC::initializeState Build Scaling <mass> type \""
169  << getDiagonalName(scaleType_) << "\"",
170  1);
171  state->invMass_ = getInvDiagonalOp(massMatrix, scaleType_);
172  }
173  // else "invMass_" should be set and there is no reason to rebuild it
174 
175  // if this is a stable discretization...we are done!
176  if (not isStabilized) {
177  state->aiD_ = Teuchos::null;
178 
179  state->setInitialized(true);
180 
181  return;
182  }
183 
184  // compute alpha scaled inv(D): EHSST2007 Eq. 4.29
185  // construct B_idF_Bt and save it for refilling later: This could reuse BQBt graph
186  LinearOp invDiagF = getInvDiagonalOp(F);
187  Teko::ModifiableLinearOp& modB_idF_Bt = state->getModifiableOp("BidFBt");
188  modB_idF_Bt = explicitMultiply(B, invDiagF, Bt, modB_idF_Bt);
189  const LinearOp B_idF_Bt = modB_idF_Bt;
190 
191  MultiVector vec_D = getDiagonal(B_idF_Bt); // this memory could be reused
192  update(-1.0, getDiagonal(C), 1.0, vec_D); // vec_D = diag(B*inv(diag(F))*Bt)-diag(C)
193  const LinearOp invD = buildInvDiagonal(vec_D, "inv(D)");
194 
195  Teko_DEBUG_MSG("Calculating alpha", 10);
196  const LinearOp BidFBtidD = multiply<double>(B_idF_Bt, invD);
197  double num = std::fabs(Teko::computeSpectralRad(BidFBtidD, 5e-2, false, eigSolveParam_));
198  Teko_DEBUG_MSG("Calculated alpha", 10);
199  state->alpha_ = 1.0 / num;
200  state->aiD_ = Thyra::scale(state->alpha_, invD);
201 
202  Teko_DEBUG_MSG_BEGIN(5) DEBUG_STREAM << "PL-LSC Alpha Parameter = " << state->alpha_ << std::endl;
203  Teko_DEBUG_MSG_END()
204 
205  state->setInitialized(true);
206 }
207 
213 void PresLaplaceLSCStrategy::computeInverses(const BlockedLinearOp& A,
214  LSCPrecondState* state) const {
215  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::computeInverses", 10);
216  Teko_DEBUG_EXPR(Teuchos::Time invTimer(""));
217 
218  std::string presLapStr = getPressureLaplaceString();
219 
220  const LinearOp F = getBlock(0, 0, A);
221  const LinearOp presLap = state->getLinearOp(presLapStr);
222 
224 
225  // (re)build the inverse of F
226  Teko_DEBUG_MSG("PL-LSC::computeInverses Building inv(F)", 1);
227  Teko_DEBUG_EXPR(invTimer.start(true));
228  ModifiableLinearOp& invF = state->getModifiableOp("invF");
229  if (invF == Teuchos::null) {
230  invF = buildInverse(*invFactoryV_, F);
231  } else {
232  rebuildInverse(*invFactoryV_, F, invF);
233  }
234  Teko_DEBUG_EXPR(invTimer.stop());
235  Teko_DEBUG_MSG("PL-LSC::computeInverses GetInvF = " << invTimer.totalElapsedTime(), 1);
236 
238 
239  // (re)build the inverse of P
240  Teko_DEBUG_MSG("PL-LSC::computeInverses Building inv(PresLap)", 1);
241  Teko_DEBUG_EXPR(invTimer.start(true));
242  ModifiableLinearOp& invPresLap = state->getModifiableOp("invPresLap");
243  if (invPresLap == Teuchos::null) {
244  invPresLap = buildInverse(*invFactoryP_, presLap);
245  } else {
246  // not need because the pressure laplacian never changes
247  // rebuildInverse(*invFactoryP_,presLap,invPresLap);
248  }
249  Teko_DEBUG_EXPR(invTimer.stop());
250  Teko_DEBUG_MSG("PL-LSC::computeInverses GetInvBQBt = " << invTimer.totalElapsedTime(), 1);
251 }
252 
254 void PresLaplaceLSCStrategy::initializeFromParameterList(const Teuchos::ParameterList& pl,
255  const InverseLibrary& invLib) {
256  // get string specifying inverse
257  std::string invStr = "", invVStr = "", invPStr = "";
258 #if defined(Teko_ENABLE_Amesos)
259  invStr = "Amesos";
260 #elif defined(Teko_ENABLE_Amesos2)
261  invStr = "Amesos2";
262 #endif
263 
264  bool useLDU = false;
265  scaleType_ = AbsRowSum;
266 
267  // "parse" the parameter list
268  if (pl.isParameter("Inverse Type")) invStr = pl.get<std::string>("Inverse Type");
269  if (pl.isParameter("Inverse Velocity Type"))
270  invVStr = pl.get<std::string>("Inverse Velocity Type");
271  if (pl.isParameter("Inverse Pressure Type"))
272  invPStr = pl.get<std::string>("Inverse Pressure Type");
273  if (pl.isParameter("Use LDU")) useLDU = pl.get<bool>("Use LDU");
274  if (pl.isParameter("Use Mass Scaling")) useMass_ = pl.get<bool>("Use Mass Scaling");
275  if (pl.isParameter("Eigen Solver Iterations"))
276  eigSolveParam_ = pl.get<int>("Eigen Solver Iterations");
277  if (pl.isParameter("Scaling Type")) {
278  scaleType_ = getDiagonalType(pl.get<std::string>("Scaling Type"));
279  TEUCHOS_TEST_FOR_EXCEPT(scaleType_ == NotDiag);
280  }
281 
282  // set defaults as needed
283  if (invVStr == "") invVStr = invStr;
284  if (invPStr == "") invPStr = invStr;
285 
286  Teko_DEBUG_MSG_BEGIN(5) DEBUG_STREAM << "LSC Inverse Strategy Parameters: " << std::endl;
287  DEBUG_STREAM << " inv v type = \"" << invVStr << "\"" << std::endl;
288  DEBUG_STREAM << " inv p type = \"" << invPStr << "\"" << std::endl;
289  DEBUG_STREAM << " use ldu = " << useLDU << std::endl;
290  DEBUG_STREAM << " use mass = " << useMass_ << std::endl;
291  DEBUG_STREAM << " scale type = " << getDiagonalName(scaleType_) << std::endl;
292  DEBUG_STREAM << "LSC Pressure Laplace Strategy Parameter list: " << std::endl;
293  pl.print(DEBUG_STREAM);
294  Teko_DEBUG_MSG_END()
295 
296  // build velocity inverse factory
297  invFactoryV_ = invLib.getInverseFactory(invVStr);
298  invFactoryP_ = invFactoryV_; // by default these are the same
299  if (invVStr != invPStr) // if different, build pressure inverse factory
300  invFactoryP_ = invLib.getInverseFactory(invPStr);
301 
302  // set other parameters
303  setUseFullLDU(useLDU);
304 }
305 
307 Teuchos::RCP<Teuchos::ParameterList> PresLaplaceLSCStrategy::getRequestedParameters() const {
308  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::getRequestedParameters", 10);
309  Teuchos::RCP<Teuchos::ParameterList> pl = rcp(new Teuchos::ParameterList());
310 
311  // grab parameters from F solver
312  RCP<Teuchos::ParameterList> fList = invFactoryV_->getRequestedParameters();
313  if (fList != Teuchos::null) {
314  Teuchos::ParameterList::ConstIterator itr;
315  for (itr = fList->begin(); itr != fList->end(); ++itr) pl->setEntry(itr->first, itr->second);
316  }
317 
318  // grab parameters from S solver
319  RCP<Teuchos::ParameterList> sList = invFactoryP_->getRequestedParameters();
320  if (sList != Teuchos::null) {
321  Teuchos::ParameterList::ConstIterator itr;
322  for (itr = sList->begin(); itr != sList->end(); ++itr) pl->setEntry(itr->first, itr->second);
323  }
324 
325  // use the mass matrix
326  if (useMass_) pl->set(getVelocityMassString(), false, "Velocity mass matrix");
327  pl->set(getPressureLaplaceString(), false, "Pressure Laplacian matrix");
328 
329  return pl;
330 }
331 
333 bool PresLaplaceLSCStrategy::updateRequestedParameters(const Teuchos::ParameterList& pl) {
334  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::updateRequestedParameters", 10);
335  bool result = true;
336 
337  // update requested parameters in solvers
338  result &= invFactoryV_->updateRequestedParameters(pl);
339  result &= invFactoryP_->updateRequestedParameters(pl);
340 
341  Teuchos::ParameterList hackList(pl);
342 
343  // get required operator acknowledgment...user must set these to true
344  bool plo = hackList.get<bool>(getPressureLaplaceString(), false);
345 
346  bool vmo = true;
347  if (useMass_) vmo = hackList.get<bool>(getVelocityMassString(), false);
348 
349  if (not plo) {
350  Teko_DEBUG_MSG("User must acknowledge the use of the \"" << getPressureLaplaceString() << "\"!",
351  0);
352  }
353  if (not vmo) {
354  Teko_DEBUG_MSG("User must acknowledge the use of the \"" << getVelocityMassString() << "\"!",
355  0);
356  }
357 
358  result &= (plo & vmo);
359 
360  return result;
361 }
362 
363 } // end namespace NS
364 } // end namespace Teko
void rebuildInverse(const InverseFactory &factory, const LinearOp &A, InverseLinearOp &invA)
An implementation of a state object for block preconditioners.
virtual Teko::LinearOp getLinearOp(const std::string &name)
Add a named operator to the state object.
LinearOp invMass_
Inverse mass operator ( )
Preconditioner state for the LSC factory.
InverseLinearOp buildInverse(const InverseFactory &factory, const LinearOp &A)
Build an inverse operator using a factory and a linear operator.
virtual void setInitialized(bool init=true)
virtual Teko::ModifiableLinearOp & getModifiableOp(const std::string &name)
Add a named operator to the state object.