Teko Version of the Day
Loading...
Searching...
No Matches
Teko_InvModALStrategy.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/*
11 * Author: Zhen Wang
12 * Email: wangz@ornl.gov
13 * zhen.wang@alum.emory.edu
14 */
15
16#include "Thyra_DefaultDiagonalLinearOp.hpp"
17#include "Thyra_DefaultAddedLinearOp_def.hpp"
18#include "Thyra_DefaultIdentityLinearOp_decl.hpp"
19
20#include "Teuchos_Time.hpp"
21
22#include "Teko_Utilities.hpp"
23
24#include "Teko_InvModALStrategy.hpp"
25#include "Teko_ModALPreconditionerFactory.hpp"
26
27using Teuchos::RCP;
28using Teuchos::rcp_const_cast;
29using Teuchos::rcp_dynamic_cast;
30
31namespace Teko {
32
33namespace NS {
34
35// Empty constructor.
36InvModALStrategy::InvModALStrategy()
37 : invFactoryA_(Teuchos::null),
38 invFactoryS_(Teuchos::null),
39 pressureMassMatrix_(Teuchos::null),
40 gamma_(0.05),
41 scaleType_(Diagonal),
42 isSymmetric_(true) {}
43
44// If there is only one InverseFactory, use it for all solves.
45InvModALStrategy::InvModALStrategy(const Teuchos::RCP<InverseFactory> &factory)
46 : invFactoryA_(factory),
47 invFactoryS_(factory),
48 pressureMassMatrix_(Teuchos::null),
49 gamma_(0.05),
50 scaleType_(Diagonal),
51 isSymmetric_(true) {}
52
53// If there are two InverseFactory's...
54InvModALStrategy::InvModALStrategy(const Teuchos::RCP<InverseFactory> &invFactoryA,
55 const Teuchos::RCP<InverseFactory> &invFactoryS)
56 : invFactoryA_(invFactoryA),
57 invFactoryS_(invFactoryS),
58 pressureMassMatrix_(Teuchos::null),
59 gamma_(0.05),
60 scaleType_(Diagonal),
61 isSymmetric_(true) {}
62
63// If there are two InverseFactory's...
64InvModALStrategy::InvModALStrategy(const Teuchos::RCP<InverseFactory> &invFactory,
65 LinearOp &pressureMassMatrix)
66 : invFactoryA_(invFactory),
67 invFactoryS_(invFactory),
68 pressureMassMatrix_(pressureMassMatrix),
69 gamma_(0.05),
70 scaleType_(Diagonal),
71 isSymmetric_(true) {}
72
73// If there are two InverseFactory's...
74InvModALStrategy::InvModALStrategy(const Teuchos::RCP<InverseFactory> &invFactoryA,
75 const Teuchos::RCP<InverseFactory> &invFactoryS,
76 LinearOp &pressureMassMatrix)
77 : invFactoryA_(invFactoryA),
78 invFactoryS_(invFactoryS),
79 pressureMassMatrix_(pressureMassMatrix),
80 gamma_(0.05),
81 scaleType_(Diagonal),
82 isSymmetric_(true) {}
83
84// Return "inverses".
85LinearOp InvModALStrategy::getInvA11p(BlockPreconditionerState &state) const {
86 return state.getInverse("invA11p");
87}
88
89LinearOp InvModALStrategy::getInvA22p(BlockPreconditionerState &state) const {
90 return state.getInverse("invA22p");
91}
92
93LinearOp InvModALStrategy::getInvA33p(BlockPreconditionerState &state) const {
94 return state.getInverse("invA33p");
95}
96
97LinearOp InvModALStrategy::getInvS(BlockPreconditionerState &state) const {
98 return state.getInverse("invS");
99}
100
101// Set pressure mass matrix.
102void InvModALStrategy::setPressureMassMatrix(const LinearOp &pressureMassMatrix) {
103 pressureMassMatrix_ = pressureMassMatrix;
104}
105
106// Set gamma.
107void InvModALStrategy::setGamma(double gamma) {
108 TEUCHOS_ASSERT(gamma > 0.0);
109 gamma_ = gamma;
110}
111
112void InvModALStrategy::buildState(const BlockedLinearOp &alOp,
113 BlockPreconditionerState &state) const {
114 Teko_DEBUG_SCOPE("InvModALStrategy::buildState", 10);
115
116 ModALPrecondState *modALState = dynamic_cast<ModALPrecondState *>(&state);
117 TEUCHOS_ASSERT(modALState != NULL);
118
119 // if necessary save state information
120 if (not modALState->isInitialized()) {
121 Teko_DEBUG_EXPR(Teuchos::Time timer(""));
122
123 {
124 // construct operators
125 Teko_DEBUG_SCOPE("ModAL::buildState: Initializing state object", 1);
126 Teko_DEBUG_EXPR(timer.start(true));
127
128 initializeState(alOp, modALState);
129
130 Teko_DEBUG_EXPR(timer.stop());
131 Teko_DEBUG_MSG("ModAL::buildState: BuildOpsTime = " << timer.totalElapsedTime(), 1);
132 }
133
134 {
135 // Build the inverses
136 Teko_DEBUG_SCOPE("ModAL::buildState: Computing inverses", 1);
137 Teko_DEBUG_EXPR(timer.start(true));
138
139 computeInverses(alOp, modALState);
140
141 Teko_DEBUG_EXPR(timer.stop());
142 Teko_DEBUG_MSG("ModAL::buildState: BuildInvTime = " << timer.totalElapsedTime(), 1);
143 }
144 }
145}
146
147// Initialize the state object using the ALOperator.
148void InvModALStrategy::initializeState(const BlockedLinearOp &alOp,
149 ModALPrecondState *state) const {
150 Teko_DEBUG_SCOPE("InvModALStrategy::initializeState", 10);
151
152 // Extract sub-matrices from blocked linear operator.
153 int dim = blockRowCount(alOp) - 1;
154 TEUCHOS_ASSERT(dim == 2 || dim == 3);
155
156 LinearOp lpA11 = getBlock(0, 0, alOp);
157 LinearOp lpA22 = getBlock(1, 1, alOp);
158 LinearOp lpA33, lpB1, lpB2, lpB3, lpB1t, lpB2t, lpB3t, lpC;
159
160 // 2D problem.
161 if (dim == 2) {
162 lpB1 = getBlock(2, 0, alOp);
163 lpB2 = getBlock(2, 1, alOp);
164 lpB1t = getBlock(0, 2, alOp);
165 lpB2t = getBlock(1, 2, alOp);
166 lpC = getBlock(2, 2, alOp);
167 }
168 // 3D problem.
169 else if (dim == 3) {
170 lpA33 = getBlock(2, 2, alOp);
171 lpB1 = getBlock(3, 0, alOp);
172 lpB2 = getBlock(3, 1, alOp);
173 lpB3 = getBlock(3, 2, alOp);
174 lpB1t = getBlock(0, 3, alOp);
175 lpB2t = getBlock(1, 3, alOp);
176 lpB3t = getBlock(2, 3, alOp);
177 lpC = getBlock(3, 3, alOp);
178 }
179
180 // For problems using stabilized finite elements,
181 // lpB1t, lpB2t and lpB3t are added linear operators. Extract original operators.
182 LinearOp B1t = (rcp_dynamic_cast<const Thyra::DefaultAddedLinearOp<double> >(lpB1t))->getOp(0);
183 LinearOp B2t = (rcp_dynamic_cast<const Thyra::DefaultAddedLinearOp<double> >(lpB2t))->getOp(0);
184 LinearOp B3t;
185 if (dim == 3) {
186 B3t = (rcp_dynamic_cast<const Thyra::DefaultAddedLinearOp<double> >(lpB3t))->getOp(0);
187 }
188
189 // std::cout << Teuchos::describe(*lpC, Teuchos::VERB_EXTREME) << std::endl;
190 // Check whether the finite elements are stable or not.
191 state->isStabilized_ = (not isZeroOp(lpC));
192 // state->isStabilized_ = false;
193 // std::cout << state->isStabilized_ << std::endl;
194
195 state->pressureMassMatrix_ = pressureMassMatrix_;
196 // If pressure mass matrix is not set, use identity.
197 if (state->pressureMassMatrix_ == Teuchos::null) {
198 Teko_DEBUG_MSG("InvModALStrategy::initializeState(): Build identity type \""
199 << getDiagonalName(scaleType_) << "\"",
200 1);
201 state->invPressureMassMatrix_ = Thyra::identity<double>(lpB1->range());
202 }
203 // If the inverse of the pressure mass matrix is not set,
204 // build it from the pressure mass matrix.
205 else if (state->invPressureMassMatrix_ == Teuchos::null) {
206 Teko_DEBUG_MSG("ModAL::initializeState(): Build Scaling <mass> type \""
207 << getDiagonalName(scaleType_) << "\"",
208 1);
209 state->invPressureMassMatrix_ = getInvDiagonalOp(pressureMassMatrix_, scaleType_);
210 }
211 // Else "invPressureMassMatrix_" should be set and there is no reason to rebuild it
212 state->gamma_ = gamma_;
213 // S_ = scale(1.0/gamma_, pressureMassMatrix_);
214 std::cout << Teuchos::describe(*(state->invPressureMassMatrix_), Teuchos::VERB_EXTREME)
215 << std::endl;
216
217 // Build state variables: B_1^T*W^{-1}*B_1, A11p, etc.
218 // B_1^T*W^{-1}*B_1 may not change so save it in the state.
219 if (state->B1tMpB1_ == Teuchos::null)
220 state->B1tMpB1_ = explicitMultiply(B1t, state->invPressureMassMatrix_, lpB1, state->B1tMpB1_);
221 // Get the(1,1) block of the non-augmented matrix.
222 // Recall alOp is augmented. So lpA11 = A11 + gamma B_1^T W^{-1} B_1.
223 // Cast lpA11 as an added linear operator and get the first item.
224 LinearOp A11 = (rcp_dynamic_cast<const Thyra::DefaultAddedLinearOp<double> >(lpA11))->getOp(0);
225 state->A11p_ = explicitAdd(A11, scale(state->gamma_, state->B1tMpB1_), state->A11p_);
226 // std::cout << Teuchos::describe(*(state->B1tMpB1_), Teuchos::VERB_EXTREME) << std::endl;
227 Teko_DEBUG_MSG("Computed A11p", 10);
228
229 if (state->B2tMpB2_ == Teuchos::null)
230 state->B2tMpB2_ = explicitMultiply(B2t, state->invPressureMassMatrix_, lpB2, state->B2tMpB2_);
231 LinearOp A22 = (rcp_dynamic_cast<const Thyra::DefaultAddedLinearOp<double> >(lpA22))->getOp(0);
232 state->A22p_ = explicitAdd(A22, scale(state->gamma_, state->B2tMpB2_), state->A22p_);
233 Teko_DEBUG_MSG("Computed A22p", 10);
234
235 if (dim == 3) {
236 if (state->B3tMpB3_ == Teuchos::null)
237 state->B3tMpB3_ = explicitMultiply(B3t, state->invPressureMassMatrix_, lpB3, state->B3tMpB3_);
238 LinearOp A33 = (rcp_dynamic_cast<const Thyra::DefaultAddedLinearOp<double> >(lpA33))->getOp(0);
239 state->A33p_ = explicitAdd(A33, scale(state->gamma_, state->B3tMpB3_), state->A33p_);
240 Teko_DEBUG_MSG("Computed A33p", 10);
241 }
242
243 // Inverse the Schur complement.
244 if (state->isStabilized_) {
245 if (state->S_ == Teuchos::null) {
246 state->S_ =
247 explicitAdd(scale(-1.0, lpC), scale(1.0 / state->gamma_, pressureMassMatrix_), state->S_);
248 }
249 Teko_DEBUG_MSG("Computed S", 10);
250 }
251
252 state->setInitialized(true);
253}
254
255// Compute inverses.
256void InvModALStrategy::computeInverses(const BlockedLinearOp &alOp,
257 ModALPrecondState *state) const {
258 int dim = blockRowCount(alOp) - 1;
259 TEUCHOS_ASSERT(dim == 2 || dim == 3);
260
261 Teko_DEBUG_SCOPE("InvModALStrategy::computeInverses", 10);
262 Teko_DEBUG_EXPR(Teuchos::Time invTimer(""));
263
264 //(re)build the inverse of A11
265 Teko_DEBUG_MSG("ModAL::computeInverses(): Building inv(A11)", 1);
266 Teko_DEBUG_EXPR(invTimer.start(true));
267
268 InverseLinearOp invA11p = state->getInverse("invA11p");
269 if (invA11p == Teuchos::null) {
270 invA11p = buildInverse(*invFactoryA_, state->A11p_);
271 state->addInverse("invA11p", invA11p);
272 } else {
273 rebuildInverse(*invFactoryA_, state->A11p_, invA11p);
274 }
275
276 Teko_DEBUG_EXPR(invTimer.stop());
277 Teko_DEBUG_MSG("ModAL::computeInverses GetInvA11 = " << invTimer.totalElapsedTime(), 1);
278
279 //(re)build the inverse of A22
280 Teko_DEBUG_MSG("ModAL::computeInverses(): Building inv(A22)", 2);
281 Teko_DEBUG_EXPR(invTimer.start(true));
282
283 InverseLinearOp invA22p = state->getInverse("invA22p");
284 if (invA22p == Teuchos::null) {
285 invA22p = buildInverse(*invFactoryA_, state->A22p_);
286 state->addInverse("invA22p", invA22p);
287 } else {
288 rebuildInverse(*invFactoryA_, state->A22p_, invA22p);
289 }
290
291 Teko_DEBUG_EXPR(invTimer.stop());
292 Teko_DEBUG_MSG("ModAL::computeInverses(): GetInvA22 = " << invTimer.totalElapsedTime(), 2);
293
294 //(re)build the inverse of A33
295 if (dim == 3) {
296 Teko_DEBUG_MSG("ModAL::computeInverses Building inv(A33)", 3);
297 Teko_DEBUG_EXPR(invTimer.start(true));
298
299 InverseLinearOp invA33p = state->getInverse("invA33p");
300 if (invA33p == Teuchos::null) {
301 invA33p = buildInverse(*invFactoryA_, state->A33p_);
302 state->addInverse("invA33p", invA33p);
303 } else {
304 rebuildInverse(*invFactoryA_, state->A33p_, invA33p);
305 }
306
307 Teko_DEBUG_EXPR(invTimer.stop());
308 Teko_DEBUG_MSG("ModAL::computeInverses GetInvA33 = " << invTimer.totalElapsedTime(), 3);
309 }
310
311 //(re)build the inverse of S
312 Teko_DEBUG_MSG("ModAL::computeInverses Building inv(S)", 4);
313 Teko_DEBUG_EXPR(invTimer.start(true));
314
315 // There are two ways to "invert" S.
316 // The following method construct invS by InverseFactory invFactoryS_.
317 // The other way is to use diagonal approximation,
318 // which is done in ModALPreconditionerFactory.cpp.
319 if (state->isStabilized_) {
320 InverseLinearOp invS = state->getInverse("invS");
321 if (invS == Teuchos::null) {
322 invS = buildInverse(*invFactoryS_, state->S_);
323 state->addInverse("invS", invS);
324 } else {
325 rebuildInverse(*invFactoryS_, state->S_, invS);
326 }
327 }
328
329 Teko_DEBUG_EXPR(invTimer.stop());
330 Teko_DEBUG_MSG("ModAL::computeInverses GetInvS = " << invTimer.totalElapsedTime(), 4);
331}
332
333} // end namespace NS
334
335} // end namespace Teko
@ Diagonal
Specifies that just the diagonal is used.
MultiVector getBlock(int i, const BlockedMultiVector &bmv)
Get the ith block from a BlockedMultiVector object.