Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_BlockInvDiagonalStrategy.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_BlockInvDiagonalStrategy.hpp"
11 #include "Teko_InverseFactory.hpp"
12 
13 namespace Teko {
14 
15 InvFactoryDiagStrategy::InvFactoryDiagStrategy(const Teuchos::RCP<InverseFactory>& factory) {
16  // only one factory to use!
17  invDiagFact_.resize(1, factory);
18  defaultInvFact_ = factory;
19 }
20 
21 InvFactoryDiagStrategy::InvFactoryDiagStrategy(
22  const std::vector<Teuchos::RCP<InverseFactory> >& factories,
23  const Teuchos::RCP<InverseFactory>& defaultFact) {
24  invDiagFact_ = factories;
25 
26  if (defaultFact == Teuchos::null)
27  defaultInvFact_ = invDiagFact_[0];
28  else
29  defaultInvFact_ = defaultFact;
30 }
31 
32 InvFactoryDiagStrategy::InvFactoryDiagStrategy(
33  const std::vector<Teuchos::RCP<InverseFactory> >& inverseFactories,
34  const std::vector<Teuchos::RCP<InverseFactory> >& preconditionerFactories,
35  const Teuchos::RCP<InverseFactory>& defaultInverseFact,
36  const Teuchos::RCP<InverseFactory>& defaultPreconditionerFact) {
37  invDiagFact_ = inverseFactories;
38  precDiagFact_ = preconditionerFactories;
39 
40  if (defaultInverseFact == Teuchos::null)
41  defaultInvFact_ = invDiagFact_[0];
42  else
43  defaultInvFact_ = defaultInverseFact;
44  defaultPrecFact_ = defaultPreconditionerFact;
45 }
46 
50 void InvFactoryDiagStrategy::getInvD(const BlockedLinearOp& A, BlockPreconditionerState& state,
51  std::vector<LinearOp>& invDiag) const {
52  Teko_DEBUG_SCOPE("InvFactoryDiagStrategy::getInvD", 10);
53 
54  // loop over diagonals, build an inverse operator for each
55  size_t diagCnt = A->productRange()->numBlocks();
56 
57  Teko_DEBUG_MSG("# diags = " << diagCnt << ", # inverses = " << invDiagFact_.size(), 6);
58 
59  const std::string opPrefix = "JacobiDiagOp";
60  for (size_t i = 0; i < diagCnt; i++) {
61  auto precFact = ((i < precDiagFact_.size()) && (!precDiagFact_[i].is_null()))
62  ? precDiagFact_[i]
63  : defaultPrecFact_;
64  auto invFact = (i < invDiagFact_.size()) ? invDiagFact_[i] : defaultInvFact_;
65  invDiag.push_back(buildInverse(*invFact, precFact, getBlock(i, i, A), state, opPrefix, i));
66  }
67 }
68 
70  Teuchos::RCP<InverseFactory>& precFact,
71  const LinearOp& matrix,
73  const std::string& opPrefix, int i) const {
74  std::stringstream ss;
75  ss << opPrefix << "_" << i;
76 
77  ModifiableLinearOp& invOp = state.getModifiableOp(ss.str());
78  ModifiableLinearOp& precOp = state.getModifiableOp("prec_" + ss.str());
79 
80  if (precFact != Teuchos::null) {
81  if (precOp == Teuchos::null) {
82  precOp = precFact->buildInverse(matrix);
83  state.addModifiableOp("prec_" + ss.str(), precOp);
84  } else {
85  Teko::rebuildInverse(*precFact, matrix, precOp);
86  }
87  }
88 
89  if (invOp == Teuchos::null)
90  if (precOp.is_null())
91  invOp = Teko::buildInverse(invFact, matrix);
92  else
93  invOp = Teko::buildInverse(invFact, matrix, precOp);
94  else {
95  if (precOp.is_null())
96  Teko::rebuildInverse(invFact, matrix, invOp);
97  else
98  Teko::rebuildInverse(invFact, matrix, precOp, invOp);
99  }
100 
101  return invOp;
102 }
103 
104 } // end namespace Teko
void rebuildInverse(const InverseFactory &factory, const LinearOp &A, InverseLinearOp &invA)
Abstract class for building an inverse operator.
virtual void getInvD(const BlockedLinearOp &A, BlockPreconditionerState &state, std::vector< LinearOp > &invDiag) const
An implementation of a state object for block preconditioners.
LinearOp buildInverse(const InverseFactory &invFact, Teuchos::RCP< InverseFactory > &precFact, const LinearOp &matrix, BlockPreconditionerState &state, const std::string &opPrefix, int i) const
Conveinence function for building inverse operators.
virtual void addModifiableOp(const std::string &name, const Teko::ModifiableLinearOp &mlo)
Add a named operator to the state object.
InverseLinearOp buildInverse(const InverseFactory &factory, const LinearOp &A)
Build an inverse operator using a factory and a linear operator.
virtual Teko::ModifiableLinearOp & getModifiableOp(const std::string &name)
Add a named operator to the state object.