Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_ModALPreconditionerFactory.cpp
1 /*
2  * Author: Zhen Wang
3  * Email: wangz@ornl.gov
4  * zhen.wang@alum.emory.edu
5  */
6 
7 #include "Teko_ModALPreconditionerFactory.hpp"
8 
9 #include "Thyra_DefaultMultipliedLinearOp.hpp"
10 #include "Thyra_DefaultAddedLinearOp.hpp"
11 #include "Thyra_DefaultIdentityLinearOp.hpp"
12 #include "Thyra_DefaultZeroLinearOp.hpp"
13 
14 #include "Teko_LU2x2InverseOp.hpp"
15 #include "Teko_Utilities.hpp"
16 #include "Teko_BlockLowerTriInverseOp.hpp"
17 #include "Teko_BlockUpperTriInverseOp.hpp"
18 #include "Teko_StaticLSCStrategy.hpp"
19 #include "Teko_InvLSCStrategy.hpp"
20 #include "Teko_PresLaplaceLSCStrategy.hpp"
21 
22 #include "Teuchos_Time.hpp"
23 
24 namespace Teko {
25 
26 namespace NS {
27 
28 ModALPrecondState::ModALPrecondState()
29  : pressureMassMatrix_(Teuchos::null), invPressureMassMatrix_(Teuchos::null) {}
30 
31 ModALPreconditionerFactory::ModALPreconditionerFactory(const Teuchos::RCP<InverseFactory>& factory)
32  : invOpsStrategy_(Teuchos::rcp(new InvModALStrategy(factory))), isSymmetric_(true) {}
33 
34 ModALPreconditionerFactory::ModALPreconditionerFactory(
35  const Teuchos::RCP<InverseFactory>& invFactoryA,
36  const Teuchos::RCP<InverseFactory>& invFactoryS)
37  : invOpsStrategy_(Teuchos::rcp(new InvModALStrategy(invFactoryA, invFactoryS))),
38  isSymmetric_(true) {}
39 
40 ModALPreconditionerFactory::ModALPreconditionerFactory(const Teuchos::RCP<InverseFactory>& factory,
41  LinearOp& pressureMassMatrix)
42  : invOpsStrategy_(Teuchos::rcp(new InvModALStrategy(factory, pressureMassMatrix))),
43  isSymmetric_(true) {}
44 
45 ModALPreconditionerFactory::ModALPreconditionerFactory(
46  const Teuchos::RCP<InverseFactory>& invFactoryA,
47  const Teuchos::RCP<InverseFactory>& invFactoryS, LinearOp& pressureMassMatrix)
48  : invOpsStrategy_(
49  Teuchos::rcp(new InvModALStrategy(invFactoryA, invFactoryS, pressureMassMatrix))),
50  isSymmetric_(true) {}
51 
52 // Construct from a strategy.
53 ModALPreconditionerFactory::ModALPreconditionerFactory(
54  const Teuchos::RCP<InvModALStrategy>& strategy)
55  : invOpsStrategy_(strategy), isSymmetric_(true) {}
56 
57 LinearOp ModALPreconditionerFactory::buildPreconditionerOperator(
58  BlockedLinearOp& alOp, BlockPreconditionerState& state) const {
59  Teko_DEBUG_SCOPE("ModALPreconditionerFactory::buildPreconditionerOperator()", 10);
60  Teko_DEBUG_EXPR(Teuchos::Time timer(""));
61  Teko_DEBUG_EXPR(Teuchos::Time totalTimer(""));
62  Teko_DEBUG_EXPR(totalTimer.start());
63 
64  // Only for 2D or 3D problems.
65  int dim = blockRowCount(alOp) - 1;
66  TEUCHOS_ASSERT(dim == 2 || dim == 3);
67 
68  // Build what is necessary for the state object.
69  Teko_DEBUG_EXPR(timer.start(true));
70  invOpsStrategy_->buildState(alOp, state);
71  Teko_DEBUG_EXPR(timer.stop());
72  Teko_DEBUG_MSG("ModALPreconditionerFactory::buildPreconditionerOperator():BuildStateTime = "
73  << timer.totalElapsedTime(),
74  2);
75 
76  // Extract inverse operators from strategy
77  Teko_DEBUG_EXPR(timer.start(true));
78  LinearOp invA11p = invOpsStrategy_->getInvA11p(state);
79  LinearOp invA22p = invOpsStrategy_->getInvA22p(state);
80  LinearOp invA33p;
81  if (dim == 3) {
82  invA33p = invOpsStrategy_->getInvA33p(state);
83  }
84 
85  // The inverse of S can be built from strategy,
86  // or just a diagonal matrix.
87  ModALPrecondState* modALState = dynamic_cast<ModALPrecondState*>(&state);
88  TEUCHOS_ASSERT(modALState != NULL);
89  LinearOp invS;
90  if (modALState->isStabilized_) {
91  invS = invOpsStrategy_->getInvS(state);
92  } else {
93  invS = scale(modALState->gamma_, modALState->invPressureMassMatrix_);
94  }
95 
96  Teko_DEBUG_EXPR(timer.stop());
97  Teko_DEBUG_MSG(
98  "ModALPrecFact::buildPreconditionerOperator(): GetInvTime = " << timer.totalElapsedTime(), 2);
99 
100  // Build diagonal operations.
101  std::vector<LinearOp> invDiag;
102  invDiag.resize(dim + 1);
103  invDiag[0] = invA11p;
104  invDiag[1] = invA22p;
105  if (dim == 2) {
106  invDiag[2] = scale(-1.0, invS);
107  } else if (dim == 3) {
108  invDiag[2] = invA33p;
109  invDiag[3] = scale(-1.0, invS);
110  }
111 
112  // Get the upper triangular matrix.
113  BlockedLinearOp U = getUpperTriBlocks(alOp);
114 
115  Teko_DEBUG_EXPR(totalTimer.stop());
116  Teko_DEBUG_MSG(
117  "ModALPrecFact::buildPreconditionerOperator TotalTime = " << totalTimer.totalElapsedTime(),
118  2);
119 
120  // Create the preconditioner
121  return createBlockUpperTriInverseOp(U, invDiag, "Modified AL preconditioner-Upper");
122 }
123 
124 } // end namespace NS
125 
126 } // end namespace Teko
An implementation of a state object for block preconditioners.
Class for saving state variables for ModALPreconditionerFactory.