Teko Version of the Day
Loading...
Searching...
No Matches
Teko_TimingsSIMPLEPreconditionerFactory.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_TimingsSIMPLEPreconditionerFactory.hpp"
11
12#include "Teko_Utilities.hpp"
13#include "Teko_InverseFactory.hpp"
14#include "Teko_BlockLowerTriInverseOp.hpp"
15#include "Teko_BlockUpperTriInverseOp.hpp"
16#ifdef TEKO_HAVE_EPETRA
17#include "Teko_DiagonalPreconditionerFactory.hpp"
18#endif
19
20#include "Teuchos_Time.hpp"
21
22using Teuchos::RCP;
23
24namespace Teko {
25namespace NS {
26
27// Constructor definition
28TimingsSIMPLEPreconditionerFactory ::TimingsSIMPLEPreconditionerFactory(
29 const RCP<InverseFactory>& inverse, double alpha)
30 : invVelFactory_(inverse),
31 invPrsFactory_(inverse),
32 alpha_(alpha),
33 fInverseType_(Diagonal),
34 useMass_(false),
35 constrTotal_("SIMPLE Constr: Total"),
36 subTotal_("SIMPLE Constr: Subs"),
37 constrCount_(0) {}
38
39TimingsSIMPLEPreconditionerFactory ::TimingsSIMPLEPreconditionerFactory(
40 const RCP<InverseFactory>& invVFact, const RCP<InverseFactory>& invPFact, double alpha)
41 : invVelFactory_(invVFact),
42 invPrsFactory_(invPFact),
43 alpha_(alpha),
44 fInverseType_(Diagonal),
45 useMass_(false),
46 constrTotal_("SIMPLE Constr: Total"),
47 subTotal_("SIMPLE Constr: Subs"),
48 constrCount_(0) {}
49
50TimingsSIMPLEPreconditionerFactory::TimingsSIMPLEPreconditionerFactory()
51 : alpha_(1.0),
52 fInverseType_(Diagonal),
53 useMass_(false),
54 constrTotal_("SIMPLE Constr: Total"),
55 subTotal_("SIMPLE Constr: Subs"),
56 constrCount_(0) {}
57
59 if (constrTotal_.totalElapsedTime() > 0.0) {
60 Teuchos::FancyOStream out(Teuchos::rcpFromRef(std::cout));
61 out.setOutputToRootOnly(0);
62
63 out << "==========================================================================="
64 << std::endl;
65 out << std::endl;
66 out << "SIMPLE Construction Count = " << constrCount_ << std::endl;
67 out << "SIMPLE Construction Total = " << constrTotal_.totalElapsedTime() << std::endl;
68 out << "SIMPLE Sub Components Total = " << subTotal_.totalElapsedTime() << std::endl;
69 out << std::endl;
70 out << "==========================================================================="
71 << std::endl;
72 }
73}
74
75// Use the factory to build the preconditioner (this is where the work goes)
76LinearOp TimingsSIMPLEPreconditionerFactory ::buildPreconditionerOperator(
77 BlockedLinearOp& blockOp, BlockPreconditionerState& state) const {
78 constrTotal_.start();
79
80 Teko_DEBUG_SCOPE("TimingsSIMPLEPreconditionerFactory::buildPreconditionerOperator", 10);
81 Teko_DEBUG_EXPR(Teuchos::Time timer(""));
82
83 int rows = blockRowCount(blockOp);
84 int cols = blockColCount(blockOp);
85
86 TEUCHOS_ASSERT(rows == 2); // sanity checks
87 TEUCHOS_ASSERT(cols == 2);
88
89 bool buildExplicitSchurComplement = true;
90
91 // extract subblocks
92 const LinearOp F = getBlock(0, 0, blockOp);
93 const LinearOp Bt = getBlock(0, 1, blockOp);
94 const LinearOp B = getBlock(1, 0, blockOp);
95 const LinearOp C = getBlock(1, 1, blockOp);
96
97 LinearOp matF = F;
98 if (useMass_) {
99 TEUCHOS_ASSERT(massMatrix_ != Teuchos::null);
100 matF = massMatrix_;
101 }
102
103 // get approximation of inv(F) name H
104 std::string fApproxStr = "<error>";
105 LinearOp H;
106 if (fInverseType_ == NotDiag) {
107 H = buildInverse(*customHFactory_, matF);
108 fApproxStr = customHFactory_->toString();
109
110 // since H is now implicit, we must build an implicit Schur complement
111 buildExplicitSchurComplement = false;
112 } else if (fInverseType_ == BlkDiag) {
113#ifdef TEKO_HAVE_EPETRA
114 // Block diagonal approximation for H
117 Hfact.initializeFromParameterList(BlkDiagList_);
118 H = Hfact.buildPreconditionerOperator(matF, Hstate);
119
120 buildExplicitSchurComplement = true; // NTS: Do I need this?
121 // Answer - no, but it is documenting whats going on here.
122#else
123 throw std::logic_error(
124 "TimingsSIMPLEPreconditionerFactory fInverseType_ == "
125 "BlkDiag but EPETRA is turned off!");
126#endif
127 } else {
128 // get generic diagonal
129 subTotal_.start();
130 H = getInvDiagonalOp(matF, fInverseType_);
131 subTotal_.stop();
132 fApproxStr = getDiagonalName(fInverseType_);
133 }
134
135 // adjust H for time scaling if it is a mass matrix
136 if (useMass_) {
137 RCP<const Teuchos::ParameterList> pl = state.getParameterList();
138
139 if (pl->isParameter("stepsize")) {
140 // get the step size
141 double stepsize = pl->get<double>("stepsize");
142
143 // scale by stepsize only if it is larger than 0
144 if (stepsize > 0.0) H = scale(stepsize, H);
145 }
146 }
147
148 // build approximate Schur complement: hatS = -C + B*H*Bt
149 LinearOp HBt, hatS;
150
151 if (buildExplicitSchurComplement) {
152 ModifiableLinearOp& mHBt = state.getModifiableOp("HBt");
153 ModifiableLinearOp& mhatS = state.getModifiableOp("hatS");
154 ModifiableLinearOp& BHBt = state.getModifiableOp("BHBt");
155
156 // build H*Bt
157 subTotal_.start();
158 mHBt = explicitMultiply(H, Bt, mHBt);
159 subTotal_.stop();
160 HBt = mHBt;
161
162 // build B*H*Bt
163 subTotal_.start();
164 BHBt = explicitMultiply(B, HBt, BHBt);
165 subTotal_.stop();
166
167 // build C-B*H*Bt
168 subTotal_.start();
169 mhatS = explicitAdd(C, scale(-1.0, BHBt), mhatS);
170 subTotal_.stop();
171 hatS = mhatS;
172 } else {
173 // build an implicit Schur complement
174 HBt = multiply(H, Bt);
175
176 hatS = add(C, scale(-1.0, multiply(B, HBt)));
177 }
178
179 // time the application of HBt
180 if (timed_HBt_ == Teuchos::null) {
181 timed_HBt_ = Teuchos::rcp(new DiagnosticLinearOp(getOutputStream(), HBt, "HBt"));
182 } else {
183 timed_HBt_->setLinearOp(HBt);
184 }
185
186 // time the application of B
187 if (timed_B_ == Teuchos::null) {
188 timed_B_ = Teuchos::rcp(new DiagnosticLinearOp(getOutputStream(), B, "B"));
189 } else {
190 timed_B_->setLinearOp(B);
191 }
192
193 // build the inverse for F
194 ModifiableLinearOp& invF = state.getModifiableOp("invF");
195 subTotal_.start();
196 if (invF == Teuchos::null) {
197 invF = buildInverse(*invVelFactory_, F);
198
199 timed_invF_ = Teuchos::rcp(new DiagnosticLinearOp(getOutputStream(), invF, "invF"));
200 } else {
201 rebuildInverse(*invVelFactory_, F, invF);
202
203 timed_invF_->setLinearOp(invF);
204 }
205 subTotal_.stop();
206
207 // build the approximate Schur complement
208 ModifiableLinearOp& invS = state.getModifiableOp("invS");
209 subTotal_.start();
210 if (invS == Teuchos::null) {
211 invS = buildInverse(*invPrsFactory_, hatS);
212
213 timed_invS_ = Teuchos::rcp(new DiagnosticLinearOp(getOutputStream(), invS, "invS"));
214 } else {
215 rebuildInverse(*invPrsFactory_, hatS, invS);
216
217 timed_invS_->setLinearOp(invS);
218 }
219 subTotal_.stop();
220
221 std::vector<LinearOp> invDiag(2); // vector storing inverses
222
223 // build lower triangular inverse matrix
224 BlockedLinearOp L = zeroBlockedOp(blockOp);
225 setBlock(1, 0, L, timed_B_);
226 endBlockFill(L);
227
228 invDiag[0] = timed_invF_;
229 invDiag[1] = timed_invS_;
230 LinearOp invL = createBlockLowerTriInverseOp(L, invDiag);
231
232 // build upper triangular matrix
233 BlockedLinearOp U = zeroBlockedOp(blockOp);
234 setBlock(0, 1, U, scale<double>(1.0 / alpha_, timed_HBt_.getConst()));
235 endBlockFill(U);
236
237 invDiag[0] = identity(rangeSpace(invF));
238 invDiag[1] = scale(alpha_, identity(rangeSpace(invS)));
239 LinearOp invU = createBlockUpperTriInverseOp(U, invDiag);
240
241 // return implicit product operator
242 Teko::LinearOp iU_t_iL = multiply(invU, invL, "SIMPLE_" + fApproxStr);
243
244 // time the application of iU_t_iL
245 if (timed_iU_t_iL_ == Teuchos::null)
246 timed_iU_t_iL_ = Teuchos::rcp(new DiagnosticLinearOp(getOutputStream(), iU_t_iL, "iU_t_iL"));
247 else
248 timed_iU_t_iL_->setLinearOp(iU_t_iL);
249
250 constrCount_++;
251
252 constrTotal_.stop();
253
254 return timed_iU_t_iL_;
255}
256
259 const Teuchos::ParameterList& pl) {
260 RCP<const InverseLibrary> invLib = getInverseLibrary();
261
262 // default conditions
263 useMass_ = false;
264 customHFactory_ = Teuchos::null;
265 fInverseType_ = Diagonal;
266
267 // get string specifying inverse
268 std::string invStr = "", invVStr = "", invPStr = "";
269 alpha_ = 1.0;
270
271 // "parse" the parameter list
272 if (pl.isParameter("Inverse Type")) invStr = pl.get<std::string>("Inverse Type");
273 if (pl.isParameter("Inverse Velocity Type"))
274 invVStr = pl.get<std::string>("Inverse Velocity Type");
275 if (pl.isParameter("Inverse Pressure Type"))
276 invPStr = pl.get<std::string>("Inverse Pressure Type");
277 if (pl.isParameter("Alpha")) alpha_ = pl.get<double>("Alpha");
278 if (pl.isParameter("Explicit Velocity Inverse Type")) {
279 std::string fInverseStr = pl.get<std::string>("Explicit Velocity Inverse Type");
280
281 // build inverse types
282 fInverseType_ = getDiagonalType(fInverseStr);
283 if (fInverseType_ == NotDiag) customHFactory_ = invLib->getInverseFactory(fInverseStr);
284
285 // Grab the sublist if we're using the block diagonal
286 if (fInverseType_ == BlkDiag) BlkDiagList_ = pl.sublist("H options");
287 }
288 if (pl.isParameter("Use Mass Scaling")) useMass_ = pl.get<bool>("Use Mass Scaling");
289
290 Teko_DEBUG_MSG_BEGIN(5) DEBUG_STREAM << "SIMPLE Parameters: " << std::endl;
291 DEBUG_STREAM << " inv type = \"" << invStr << "\"" << std::endl;
292 DEBUG_STREAM << " inv v type = \"" << invVStr << "\"" << std::endl;
293 DEBUG_STREAM << " inv p type = \"" << invPStr << "\"" << std::endl;
294 DEBUG_STREAM << " alpha = " << alpha_ << std::endl;
295 DEBUG_STREAM << " use mass = " << useMass_ << std::endl;
296 DEBUG_STREAM << " vel scaling = " << getDiagonalName(fInverseType_) << std::endl;
297 DEBUG_STREAM << "SIMPLE Parameter list: " << std::endl;
298 pl.print(DEBUG_STREAM);
299 Teko_DEBUG_MSG_END()
300
301 // set defaults as needed
302#if defined(Teko_ENABLE_Amesos)
303 if (invStr == "") invStr = "Amesos";
304#elif defined(Teko_ENABLE_Amesos2)
305 if (invStr == "") invStr = "Amesos2";
306#endif
307 if (invVStr == "") invVStr = invStr;
308 if (invPStr == "") invPStr = invStr;
309
310 // two inverse factory objects
311 RCP<InverseFactory> invVFact, invPFact;
312
313 // build velocity inverse factory
314 invVFact = invLib->getInverseFactory(invVStr);
315 invPFact = invVFact; // by default these are the same
316 if (invVStr != invPStr) // if different, build pressure inverse factory
317 invPFact = invLib->getInverseFactory(invPStr);
318
319 // based on parameter type build a strategy
320 invVelFactory_ = invVFact;
321 invPrsFactory_ = invPFact;
322
323 if (useMass_) {
324 Teuchos::RCP<Teko::RequestHandler> rh = getRequestHandler();
325 rh->preRequest<Teko::LinearOp>(Teko::RequestMesg("Velocity Mass Matrix"));
326 Teko::LinearOp mass = rh->request<Teko::LinearOp>(Teko::RequestMesg("Velocity Mass Matrix"));
327 setMassMatrix(mass);
328 }
329}
330
333 const {
334 Teuchos::RCP<Teuchos::ParameterList> result;
335 Teuchos::RCP<Teuchos::ParameterList> pl = rcp(new Teuchos::ParameterList());
336
337 // grab parameters from F solver
338 RCP<Teuchos::ParameterList> vList = invVelFactory_->getRequestedParameters();
339 if (vList != Teuchos::null) {
340 Teuchos::ParameterList::ConstIterator itr;
341 for (itr = vList->begin(); itr != vList->end(); ++itr) pl->setEntry(itr->first, itr->second);
342 result = pl;
343 }
344
345 // grab parameters from S solver
346 RCP<Teuchos::ParameterList> pList = invPrsFactory_->getRequestedParameters();
347 if (pList != Teuchos::null) {
348 Teuchos::ParameterList::ConstIterator itr;
349 for (itr = pList->begin(); itr != pList->end(); ++itr) pl->setEntry(itr->first, itr->second);
350 result = pl;
351 }
352
353 // grab parameters from S solver
354 if (customHFactory_ != Teuchos::null) {
355 RCP<Teuchos::ParameterList> hList = customHFactory_->getRequestedParameters();
356 if (hList != Teuchos::null) {
357 Teuchos::ParameterList::ConstIterator itr;
358 for (itr = hList->begin(); itr != hList->end(); ++itr) pl->setEntry(itr->first, itr->second);
359 result = pl;
360 }
361 }
362
363 return result;
364}
365
368 const Teuchos::ParameterList& pl) {
369 Teko_DEBUG_SCOPE("InvLSCStrategy::updateRequestedParameters", 10);
370 bool result = true;
371
372 // update requested parameters in solvers
373 result &= invVelFactory_->updateRequestedParameters(pl);
374 result &= invPrsFactory_->updateRequestedParameters(pl);
375 if (customHFactory_ != Teuchos::null) result &= customHFactory_->updateRequestedParameters(pl);
376
377 return result;
378}
379
380} // end namespace NS
381} // end namespace Teko
@ BlkDiag
Specifies that a block diagonal approximation is to be used.
@ 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.
This linear operator prints diagnostics about operator application and creation times....
Preconditioner factory for building explcit inverse of diagonal operators. This includes block operat...
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl)
Initialize from a parameter list.
LinearOp buildPreconditionerOperator(LinearOp &lo, PreconditionerState &state) const
virtual bool updateRequestedParameters(const Teuchos::ParameterList &pl)
For assisting in construction of the preconditioner.
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl)
Initialize from a parameter list.
virtual void setMassMatrix(Teko::LinearOp &mass)
Set the mass matrix for this factory.
virtual ~TimingsSIMPLEPreconditionerFactory()
Destructor that outputs construction timings.
virtual Teuchos::RCP< Teuchos::ParameterList > getRequestedParameters() const
For assisting in construction of the preconditioner.
Teuchos::RCP< const InverseLibrary > getInverseLibrary() const
Get the inverse library used by this preconditioner factory.
Teuchos::RCP< RequestHandler > getRequestHandler() const
Get the request handler with pointers to the appropriate callbacks.
virtual Teko::ModifiableLinearOp & getModifiableOp(const std::string &name)
Add a named operator to the state object.