Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_PresLaplaceLSCStrategy.cpp
1 /*
2 // @HEADER
3 //
4 // ***********************************************************************
5 //
6 // Teko: A package for block and physics based preconditioning
7 // Copyright 2010 Sandia Corporation
8 //
9 // Under the terms of Contract DE-AC04-94AL85000 with Sandia Corporation,
10 // the U.S. Government retains certain rights in this software.
11 //
12 // Redistribution and use in source and binary forms, with or without
13 // modification, are permitted provided that the following conditions are
14 // met:
15 //
16 // 1. Redistributions of source code must retain the above copyright
17 // notice, this list of conditions and the following disclaimer.
18 //
19 // 2. Redistributions in binary form must reproduce the above copyright
20 // notice, this list of conditions and the following disclaimer in the
21 // documentation and/or other materials provided with the distribution.
22 //
23 // 3. Neither the name of the Corporation nor the names of the
24 // contributors may be used to endorse or promote products derived from
25 // this software without specific prior written permission.
26 //
27 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
28 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
29 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
30 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
31 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
32 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
33 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
34 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
35 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
36 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
37 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
38 //
39 // Questions? Contact Eric C. Cyr (eccyr@sandia.gov)
40 //
41 // ***********************************************************************
42 //
43 // @HEADER
44 
45 */
46 
47 #include "NS/Teko_PresLaplaceLSCStrategy.hpp"
48 
49 #include "Thyra_DefaultDiagonalLinearOp.hpp"
50 #include "Thyra_EpetraThyraWrappers.hpp"
51 #include "Thyra_get_Epetra_Operator.hpp"
52 #include "Thyra_EpetraLinearOp.hpp"
53 
54 #include "Teuchos_Time.hpp"
55 #include "Teuchos_TimeMonitor.hpp"
56 
57 // Teko includes
58 #include "Teko_Utilities.hpp"
59 #include "NS/Teko_LSCPreconditionerFactory.hpp"
60 
61 using Teuchos::RCP;
62 using Teuchos::rcp_dynamic_cast;
63 using Teuchos::rcp_const_cast;
64 
65 namespace Teko {
66 namespace NS {
67 
69 // PresLaplaceLSCStrategy Implementation
71 
72 // constructors
74 PresLaplaceLSCStrategy::PresLaplaceLSCStrategy()
75  : invFactoryV_(Teuchos::null), invFactoryP_(Teuchos::null)
76  , eigSolveParam_(5), useFullLDU_(false), useMass_(false), scaleType_(AbsRowSum)
77 { }
78 
79 PresLaplaceLSCStrategy::PresLaplaceLSCStrategy(const Teuchos::RCP<InverseFactory> & factory)
80  : invFactoryV_(factory), invFactoryP_(factory)
81  , eigSolveParam_(5), useFullLDU_(false), useMass_(false), scaleType_(AbsRowSum)
82 { }
83 
84 PresLaplaceLSCStrategy::PresLaplaceLSCStrategy(const Teuchos::RCP<InverseFactory> & invFactF,
85  const Teuchos::RCP<InverseFactory> & invFactS)
86  : invFactoryV_(invFactF), invFactoryP_(invFactS)
87  , eigSolveParam_(5), useFullLDU_(false), useMass_(false), scaleType_(AbsRowSum)
88 { }
89 
91 
92 void PresLaplaceLSCStrategy::buildState(BlockedLinearOp & A,BlockPreconditionerState & state) const
93 {
94  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::buildState",10);
95 
96  LSCPrecondState * lscState = dynamic_cast<LSCPrecondState*>(&state);
97  TEUCHOS_ASSERT(lscState!=0);
98 
99  // if neccessary save state information
100  if(not lscState->isInitialized()) {
101  Teko_DEBUG_EXPR(Teuchos::Time timer(""));
102 
103  // construct operators
104  {
105  Teko_DEBUG_SCOPE("PL-LSC::buildState constructing operators",1);
106  Teko_DEBUG_EXPR(timer.start(true));
107 
108  initializeState(A,lscState);
109 
110  Teko_DEBUG_EXPR(timer.stop());
111  Teko_DEBUG_MSG("PL-LSC::buildState BuildOpsTime = " << timer.totalElapsedTime(),1);
112  }
113 
114  // Build the inverses
115  {
116  Teko_DEBUG_SCOPE("PL-LSC::buildState calculating inverses",1);
117  Teko_DEBUG_EXPR(timer.start(true));
118 
119  computeInverses(A,lscState);
120 
121  Teko_DEBUG_EXPR(timer.stop());
122  Teko_DEBUG_MSG("PL-LSC::buildState BuildInvTime = " << timer.totalElapsedTime(),1);
123  }
124  }
125 }
126 
127 // functions inherited from LSCStrategy
128 LinearOp PresLaplaceLSCStrategy::getInvBQBt(const BlockedLinearOp & /* A */,BlockPreconditionerState & state) const
129 {
130  return state.getModifiableOp("invPresLap");
131 }
132 
133 LinearOp PresLaplaceLSCStrategy::getInvBHBt(const BlockedLinearOp & /* A */,BlockPreconditionerState & state) const
134 {
135  return state.getModifiableOp("invPresLap");
136 }
137 
138 LinearOp PresLaplaceLSCStrategy::getInvF(const BlockedLinearOp & /* A */,BlockPreconditionerState & state) const
139 {
140  return state.getModifiableOp("invF");
141 }
142 
143 // LinearOp PresLaplaceLSCStrategy::getInvAlphaD(const BlockedLinearOp & A,BlockPreconditionerState & state) const
144 LinearOp PresLaplaceLSCStrategy::getOuterStabilization(const BlockedLinearOp & /* A */,BlockPreconditionerState & state) const
145 {
146  LSCPrecondState * lscState = dynamic_cast<LSCPrecondState*>(&state);
147  TEUCHOS_ASSERT(lscState!=0);
148  TEUCHOS_ASSERT(lscState->isInitialized())
149 
150  return lscState->aiD_;
151 }
152 
153 LinearOp PresLaplaceLSCStrategy::getInvMass(const BlockedLinearOp & /* A */,BlockPreconditionerState & state) const
154 {
155  LSCPrecondState * lscState = dynamic_cast<LSCPrecondState*>(&state);
156  TEUCHOS_ASSERT(lscState!=0);
157  TEUCHOS_ASSERT(lscState->isInitialized())
158 
159  return lscState->invMass_;
160 }
161 
162 LinearOp PresLaplaceLSCStrategy::getHScaling(const BlockedLinearOp & A,BlockPreconditionerState & state) const
163 {
164  return getInvMass(A,state);
165 }
166 
168 void PresLaplaceLSCStrategy::initializeState(const BlockedLinearOp & A,LSCPrecondState * state) const
169 {
170  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::initializeState",10);
171 
172  std::string velMassStr = getVelocityMassString();
173 
174  const LinearOp F = getBlock(0,0,A);
175  const LinearOp Bt = getBlock(0,1,A);
176  const LinearOp B = getBlock(1,0,A);
177  const LinearOp C = getBlock(1,1,A);
178 
179  LinearOp D = B;
180  LinearOp G = Bt;
181 
182  bool isStabilized = (not isZeroOp(C));
183 
184  // grab operators from state object
185  LinearOp massMatrix = state->getLinearOp(velMassStr);
186 
187  // The logic follows like this
188  // if there is no mass matrix available --> build from F
189  // if there is a mass matrix and the inverse hasn't yet been built
190  // --> build from the mass matrix
191  // otherwise, there is already an invMass_ matrix that is appropriate
192  // --> use that one
193  if(massMatrix==Teuchos::null) {
194  Teko_DEBUG_MSG("PL-LSC::initializeState Build Scaling <F> type \""
195  << getDiagonalName(scaleType_) << "\"" ,1);
196  state->invMass_ = getInvDiagonalOp(F,scaleType_);
197  }
198  else if(state->invMass_==Teuchos::null) {
199  Teko_DEBUG_MSG("PL-LSC::initializeState Build Scaling <mass> type \""
200  << getDiagonalName(scaleType_) << "\"" ,1);
201  state->invMass_ = getInvDiagonalOp(massMatrix,scaleType_);
202  }
203  // else "invMass_" should be set and there is no reason to rebuild it
204 
205  // if this is a stable discretization...we are done!
206  if(not isStabilized) {
207  state->aiD_ = Teuchos::null;
208 
209  state->setInitialized(true);
210 
211  return;
212  }
213 
214  // compute alpha scaled inv(D): EHSST2007 Eq. 4.29
215  // construct B_idF_Bt and save it for refilling later: This could reuse BQBt graph
216  LinearOp invDiagF = getInvDiagonalOp(F);
217  Teko::ModifiableLinearOp & modB_idF_Bt = state->getModifiableOp("BidFBt");
218  modB_idF_Bt = explicitMultiply(B,invDiagF,Bt,modB_idF_Bt);
219  const LinearOp B_idF_Bt = modB_idF_Bt;
220 
221  MultiVector vec_D = getDiagonal(B_idF_Bt); // this memory could be reused
222  update(-1.0,getDiagonal(C),1.0,vec_D); // vec_D = diag(B*inv(diag(F))*Bt)-diag(C)
223  const LinearOp invD = buildInvDiagonal(vec_D,"inv(D)");
224 
225  Teko_DEBUG_MSG("Calculating alpha",10);
226  const LinearOp BidFBtidD = multiply<double>(B_idF_Bt,invD);
227  double num = std::fabs(Teko::computeSpectralRad(BidFBtidD,5e-2,false,eigSolveParam_));
228  Teko_DEBUG_MSG("Calculated alpha",10);
229  state->alpha_ = 1.0/num;
230  state->aiD_ = Thyra::scale(state->alpha_,invD);
231 
232  Teko_DEBUG_MSG_BEGIN(5)
233  DEBUG_STREAM << "PL-LSC Alpha Parameter = " << state->alpha_ << std::endl;
234  Teko_DEBUG_MSG_END()
235 
236  state->setInitialized(true);
237 }
238 
244 void PresLaplaceLSCStrategy::computeInverses(const BlockedLinearOp & A,LSCPrecondState * state) const
245 {
246  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::computeInverses",10);
247  Teko_DEBUG_EXPR(Teuchos::Time invTimer(""));
248 
249  std::string presLapStr = getPressureLaplaceString();
250 
251  const LinearOp F = getBlock(0,0,A);
252  const LinearOp presLap = state->getLinearOp(presLapStr);
253 
255 
256  // (re)build the inverse of F
257  Teko_DEBUG_MSG("PL-LSC::computeInverses Building inv(F)",1);
258  Teko_DEBUG_EXPR(invTimer.start(true));
259  ModifiableLinearOp & invF = state->getModifiableOp("invF");
260  if(invF==Teuchos::null) {
261  invF = buildInverse(*invFactoryV_,F);
262  } else {
263  rebuildInverse(*invFactoryV_,F,invF);
264  }
265  Teko_DEBUG_EXPR(invTimer.stop());
266  Teko_DEBUG_MSG("PL-LSC::computeInverses GetInvF = " << invTimer.totalElapsedTime(),1);
267 
269 
270  // (re)build the inverse of P
271  Teko_DEBUG_MSG("PL-LSC::computeInverses Building inv(PresLap)",1);
272  Teko_DEBUG_EXPR(invTimer.start(true));
273  ModifiableLinearOp & invPresLap = state->getModifiableOp("invPresLap");
274  if(invPresLap==Teuchos::null) {
275  invPresLap = buildInverse(*invFactoryP_,presLap);
276  } else {
277  // not need because the pressure laplacian never changes
278  // rebuildInverse(*invFactoryP_,presLap,invPresLap);
279  }
280  Teko_DEBUG_EXPR(invTimer.stop());
281  Teko_DEBUG_MSG("PL-LSC::computeInverses GetInvBQBt = " << invTimer.totalElapsedTime(),1);
282 }
283 
285 void PresLaplaceLSCStrategy::initializeFromParameterList(const Teuchos::ParameterList & pl,const InverseLibrary & invLib)
286 {
287  // get string specifying inverse
288  std::string invStr="Amesos", invVStr="", invPStr="";
289  bool useLDU = false;
290  scaleType_ = AbsRowSum;
291 
292  // "parse" the parameter list
293  if(pl.isParameter("Inverse Type"))
294  invStr = pl.get<std::string>("Inverse Type");
295  if(pl.isParameter("Inverse Velocity Type"))
296  invVStr = pl.get<std::string>("Inverse Velocity Type");
297  if(pl.isParameter("Inverse Pressure Type"))
298  invPStr = pl.get<std::string>("Inverse Pressure Type");
299  if(pl.isParameter("Use LDU"))
300  useLDU = pl.get<bool>("Use LDU");
301  if(pl.isParameter("Use Mass Scaling"))
302  useMass_ = pl.get<bool>("Use Mass Scaling");
303  if(pl.isParameter("Eigen Solver Iterations"))
304  eigSolveParam_ = pl.get<int>("Eigen Solver Iterations");
305  if(pl.isParameter("Scaling Type")) {
306  scaleType_ = getDiagonalType(pl.get<std::string>("Scaling Type"));
307  TEUCHOS_TEST_FOR_EXCEPT(scaleType_==NotDiag);
308  }
309 
310  // set defaults as needed
311  if(invVStr=="") invVStr = invStr;
312  if(invPStr=="") invPStr = invStr;
313 
314  Teko_DEBUG_MSG_BEGIN(5)
315  DEBUG_STREAM << "LSC Inverse Strategy Parameters: " << std::endl;
316  DEBUG_STREAM << " inv v type = \"" << invVStr << "\"" << std::endl;
317  DEBUG_STREAM << " inv p type = \"" << invPStr << "\"" << std::endl;
318  DEBUG_STREAM << " use ldu = " << useLDU << std::endl;
319  DEBUG_STREAM << " use mass = " << useMass_ << std::endl;
320  DEBUG_STREAM << " scale type = " << getDiagonalName(scaleType_) << std::endl;
321  DEBUG_STREAM << "LSC Pressure Laplace Strategy Parameter list: " << std::endl;
322  pl.print(DEBUG_STREAM);
323  Teko_DEBUG_MSG_END()
324 
325  // build velocity inverse factory
326  invFactoryV_ = invLib.getInverseFactory(invVStr);
327  invFactoryP_ = invFactoryV_; // by default these are the same
328  if(invVStr!=invPStr) // if different, build pressure inverse factory
329  invFactoryP_ = invLib.getInverseFactory(invPStr);
330 
331  // set other parameters
332  setUseFullLDU(useLDU);
333 }
334 
336 Teuchos::RCP<Teuchos::ParameterList> PresLaplaceLSCStrategy::getRequestedParameters() const
337 {
338  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::getRequestedParameters",10);
339  Teuchos::RCP<Teuchos::ParameterList> pl = rcp(new Teuchos::ParameterList());
340 
341  // grab parameters from F solver
342  RCP<Teuchos::ParameterList> fList = invFactoryV_->getRequestedParameters();
343  if(fList!=Teuchos::null) {
344  Teuchos::ParameterList::ConstIterator itr;
345  for(itr=fList->begin();itr!=fList->end();++itr)
346  pl->setEntry(itr->first,itr->second);
347  }
348 
349  // grab parameters from S solver
350  RCP<Teuchos::ParameterList> sList = invFactoryP_->getRequestedParameters();
351  if(sList!=Teuchos::null) {
352  Teuchos::ParameterList::ConstIterator itr;
353  for(itr=sList->begin();itr!=sList->end();++itr)
354  pl->setEntry(itr->first,itr->second);
355  }
356 
357  // use the mass matrix
358  if(useMass_)
359  pl->set(getVelocityMassString(), false,"Velocity mass matrix");
360  pl->set(getPressureLaplaceString(), false,"Pressure Laplacian matrix");
361 
362  return pl;
363 }
364 
366 bool PresLaplaceLSCStrategy::updateRequestedParameters(const Teuchos::ParameterList & pl)
367 {
368  Teko_DEBUG_SCOPE("PresLaplaceLSCStrategy::updateRequestedParameters",10);
369  bool result = true;
370 
371  // update requested parameters in solvers
372  result &= invFactoryV_->updateRequestedParameters(pl);
373  result &= invFactoryP_->updateRequestedParameters(pl);
374 
375  Teuchos::ParameterList hackList(pl);
376 
377  // get required operator acknowledgment...user must set these to true
378  bool plo = hackList.get<bool>(getPressureLaplaceString(),false);
379 
380  bool vmo = true;
381  if(useMass_)
382  vmo = hackList.get<bool>(getVelocityMassString(),false);
383 
384  if(not plo) { Teko_DEBUG_MSG("User must acknowledge the use of the \"" << getPressureLaplaceString() << "\"!",0); }
385  if(not vmo) { Teko_DEBUG_MSG("User must acknowledge the use of the \"" << getVelocityMassString() << "\"!",0); }
386 
387  result &= (plo & vmo);
388 
389  return result;
390 }
391 
392 } // end namespace NS
393 } // 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.