Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_LSCSIMPLECStrategy.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 "Teko_LSCSIMPLECStrategy.hpp"
48 
49 #include "Thyra_DefaultDiagonalLinearOp.hpp"
50 #include "Thyra_VectorStdOps.hpp"
51 
52 #include "Teuchos_Time.hpp"
53 #include "Teuchos_TimeMonitor.hpp"
54 
55 // Teko includes
56 #include "Teko_Utilities.hpp"
57 #include "Teko_LSCPreconditionerFactory.hpp"
58 
59 using Teuchos::RCP;
60 using Teuchos::rcp_const_cast;
61 using Teuchos::rcp_dynamic_cast;
62 
63 namespace Teko {
64 namespace NS {
65 
67 // LSCSIMPLECStrategy Implementation
69 
70 // constructors
72 LSCSIMPLECStrategy::LSCSIMPLECStrategy()
73  : invFactoryF_(Teuchos::null),
74  invFactoryS_(Teuchos::null),
75  useFullLDU_(false),
76  scaleType_(Diagonal) {}
77 
79 
80 void LSCSIMPLECStrategy::buildState(BlockedLinearOp& A, BlockPreconditionerState& state) const {
81  Teko_DEBUG_SCOPE("LSCSIMPLECStrategy::buildState", 10);
82 
83  LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
84  TEUCHOS_ASSERT(lscState != 0);
85 
86  // if neccessary save state information
87  if (not lscState->isInitialized()) {
88  Teko_DEBUG_EXPR(Teuchos::Time timer(""));
89 
90  // construct operators
91  {
92  Teko_DEBUG_SCOPE("LSC-SIMPLEC::buildState constructing operators", 1);
93  Teko_DEBUG_EXPR(timer.start(true));
94 
95  initializeState(A, lscState);
96 
97  Teko_DEBUG_EXPR(timer.stop());
98  Teko_DEBUG_MSG("LSC-SIMPLEC::buildState BuildOpsTime = " << timer.totalElapsedTime(), 1);
99  }
100 
101  // Build the inverses
102  {
103  Teko_DEBUG_SCOPE("LSC-SIMPLEC::buildState calculating inverses", 1);
104  Teko_DEBUG_EXPR(timer.start(true));
105 
106  computeInverses(A, lscState);
107 
108  Teko_DEBUG_EXPR(timer.stop());
109  Teko_DEBUG_MSG("LSC-SIMPLEC::buildState BuildInvTime = " << timer.totalElapsedTime(), 1);
110  }
111  }
112 }
113 
114 // functions inherited from LSCStrategy
115 LinearOp LSCSIMPLECStrategy::getInvBQBt(const BlockedLinearOp& /* A */,
116  BlockPreconditionerState& state) const {
117  return state.getInverse("invBQBtmC");
118 }
119 
120 LinearOp LSCSIMPLECStrategy::getInvBHBt(const BlockedLinearOp& /* A */,
121  BlockPreconditionerState& state) const {
122  return state.getInverse("invBQBtmC").getConst();
123 }
124 
125 LinearOp LSCSIMPLECStrategy::getInvF(const BlockedLinearOp& /* A */,
126  BlockPreconditionerState& state) const {
127  return state.getInverse("invF");
128 }
129 
130 LinearOp LSCSIMPLECStrategy::getInnerStabilization(const BlockedLinearOp& A,
131  BlockPreconditionerState& state) const {
132  LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
133  TEUCHOS_ASSERT(lscState != 0);
134  TEUCHOS_ASSERT(lscState->isInitialized())
135 
136  const LinearOp C = getBlock(1, 1, A);
137  return scale(-1.0, C);
138 }
139 
140 LinearOp LSCSIMPLECStrategy::getInvMass(const BlockedLinearOp& /* A */,
141  BlockPreconditionerState& state) const {
142  LSCPrecondState* lscState = dynamic_cast<LSCPrecondState*>(&state);
143  TEUCHOS_ASSERT(lscState != 0);
144  TEUCHOS_ASSERT(lscState->isInitialized())
145 
146  return lscState->invMass_;
147 }
148 
149 LinearOp LSCSIMPLECStrategy::getHScaling(const BlockedLinearOp& A,
150  BlockPreconditionerState& state) const {
151  return getInvMass(A, state);
152 }
153 
155 void LSCSIMPLECStrategy::initializeState(const BlockedLinearOp& A, LSCPrecondState* state) const {
156  Teko_DEBUG_SCOPE("LSCSIMPLECStrategy::initializeState", 10);
157 
158  const LinearOp F = getBlock(0, 0, A);
159  const LinearOp Bt = getBlock(0, 1, A);
160  const LinearOp B = getBlock(1, 0, A);
161  const LinearOp C = getBlock(1, 1, A);
162 
163  bool isStabilized = (not isZeroOp(C));
164 
165  state->invMass_ = getInvDiagonalOp(F, scaleType_);
166 
167  // compute BQBt
168  state->BQBt_ = explicitMultiply(B, state->invMass_, Bt, state->BQBt_);
169  if (isStabilized) {
170  // now build B*Q*Bt-C
171  Teko::ModifiableLinearOp BQBtmC = state->getInverse("BQBtmC");
172  BQBtmC = explicitAdd(state->BQBt_, scale(-1.0, C), BQBtmC);
173  state->addInverse("BQBtmC", BQBtmC);
174  }
175  Teko_DEBUG_MSG("Computed BQBt", 10);
176 
177  state->setInitialized(true);
178 }
179 
185 void LSCSIMPLECStrategy::computeInverses(const BlockedLinearOp& A, LSCPrecondState* state) const {
186  Teko_DEBUG_SCOPE("LSCSIMPLECStrategy::computeInverses", 10);
187  Teko_DEBUG_EXPR(Teuchos::Time invTimer(""));
188 
189  const LinearOp F = getBlock(0, 0, A);
190 
192 
193  // (re)build the inverse of F
194  Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses Building inv(F)", 1);
195  Teko_DEBUG_EXPR(invTimer.start(true));
196  InverseLinearOp invF = state->getInverse("invF");
197  if (invF == Teuchos::null) {
198  invF = buildInverse(*invFactoryF_, F);
199  state->addInverse("invF", invF);
200  } else {
201  rebuildInverse(*invFactoryF_, F, invF);
202  }
203  Teko_DEBUG_EXPR(invTimer.stop());
204  Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses GetInvF = " << invTimer.totalElapsedTime(), 1);
205 
207 
208  // (re)build the inverse of BQBt
209  Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses Building inv(BQBtmC)", 1);
210  Teko_DEBUG_EXPR(invTimer.start(true));
211  const LinearOp BQBt = state->getInverse("BQBtmC");
212  InverseLinearOp invBQBt = state->getInverse("invBQBtmC");
213  if (invBQBt == Teuchos::null) {
214  invBQBt = buildInverse(*invFactoryS_, BQBt);
215  state->addInverse("invBQBtmC", invBQBt);
216  } else {
217  rebuildInverse(*invFactoryS_, BQBt, invBQBt);
218  }
219  Teko_DEBUG_EXPR(invTimer.stop());
220  Teko_DEBUG_MSG("LSC-SIMPLEC::computeInverses GetInvBQBt = " << invTimer.totalElapsedTime(), 1);
221 }
222 
224 void LSCSIMPLECStrategy::initializeFromParameterList(const Teuchos::ParameterList& pl,
225  const InverseLibrary& invLib) {
226  // get string specifying inverse
227  std::string invStr = "", invVStr = "", invPStr = "";
228  bool useLDU = false;
229  scaleType_ = Diagonal;
230 
231  // "parse" the parameter list
232  if (pl.isParameter("Inverse Type")) invStr = pl.get<std::string>("Inverse Type");
233  if (pl.isParameter("Inverse Velocity Type"))
234  invVStr = pl.get<std::string>("Inverse Velocity Type");
235  if (pl.isParameter("Inverse Pressure Type"))
236  invPStr = pl.get<std::string>("Inverse Pressure Type");
237  if (pl.isParameter("Use LDU")) useLDU = pl.get<bool>("Use LDU");
238  if (pl.isParameter("Scaling Type")) {
239  scaleType_ = getDiagonalType(pl.get<std::string>("Scaling Type"));
240  TEUCHOS_TEST_FOR_EXCEPT(scaleType_ == NotDiag);
241  }
242 
243  Teko_DEBUG_MSG_BEGIN(0) DEBUG_STREAM << "LSC Inverse Strategy Parameters: " << std::endl;
244  DEBUG_STREAM << " inv type = \"" << invStr << "\"" << std::endl;
245  DEBUG_STREAM << " inv v type = \"" << invVStr << "\"" << std::endl;
246  DEBUG_STREAM << " inv p type = \"" << invPStr << "\"" << std::endl;
247  DEBUG_STREAM << " use ldu = " << useLDU << std::endl;
248  DEBUG_STREAM << " scale type = " << getDiagonalName(scaleType_) << std::endl;
249  DEBUG_STREAM << "LSC Inverse Strategy Parameter list: " << std::endl;
250  pl.print(DEBUG_STREAM);
251  Teko_DEBUG_MSG_END()
252 
253  // set defaults as needed
254 #if defined(Teko_ENABLE_Amesos)
255  if (invStr == "") invStr = "Amesos";
256 #elif defined(Teko_ENABLE_Amesos2)
257  if (invStr == "") invStr = "Amesos2";
258 #endif
259  if (invVStr == "") invVStr = invStr;
260  if (invPStr == "") invPStr = invStr;
261 
262  // build velocity inverse factory
263  invFactoryF_ = invLib.getInverseFactory(invVStr);
264  invFactoryS_ = invFactoryF_; // by default these are the same
265  if (invVStr != invPStr) // if different, build pressure inverse factory
266  invFactoryS_ = invLib.getInverseFactory(invPStr);
267 
268  // set other parameters
269  setUseFullLDU(useLDU);
270 }
271 
272 } // end namespace NS
273 } // end namespace Teko
void rebuildInverse(const InverseFactory &factory, const LinearOp &A, InverseLinearOp &invA)
virtual Teko::InverseLinearOp getInverse(const std::string &name) const
Get a named inverse from the state object.
void computeInverses(const BlockedLinearOp &A, LSCPrecondState *state) const
virtual LinearOp getInvMass(const BlockedLinearOp &A, BlockPreconditionerState &state) const
An implementation of a state object for block preconditioners.
virtual void setUseFullLDU(bool val)
Set to true to use the Full LDU decomposition, false otherwise.
virtual void addInverse(const std::string &name, const Teko::InverseLinearOp &ilo)
Add a named inverse to the state object.
virtual void initializeState(const BlockedLinearOp &A, LSCPrecondState *state) const
Initialize the state object using this blocked linear operator.
virtual LinearOp getHScaling(const BlockedLinearOp &A, BlockPreconditionerState &state) const
LinearOp invMass_
Inverse mass operator ( )
Preconditioner state for the LSC factory.
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl, const InverseLibrary &invLib)
Initialize from a parameter list.
virtual LinearOp getInnerStabilization(const BlockedLinearOp &A, BlockPreconditionerState &state) const
InverseLinearOp buildInverse(const InverseFactory &factory, const LinearOp &A)
Build an inverse operator using a factory and a linear operator.
virtual LinearOp getInvF(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual LinearOp getInvBHBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual LinearOp getInvBQBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void buildState(BlockedLinearOp &A, BlockPreconditionerState &state) const
Functions inherited from LSCStrategy.
virtual void setInitialized(bool init=true)