Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_JacobiPreconditionerFactory.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_JacobiPreconditionerFactory.hpp"
11 #include "Teko_BlockDiagonalInverseOp.hpp"
12 
13 using Teuchos::rcp;
14 
15 namespace Teko {
16 
18  const LinearOp& invD1)
19  : invOpsStrategy_(rcp(new StaticInvDiagStrategy(invD0, invD1))) {}
20 
22  const RCP<const BlockInvDiagonalStrategy>& strategy)
23  : invOpsStrategy_(strategy) {}
24 
28 
30  BlockedLinearOp& blo, BlockPreconditionerState& state) const {
31  int rows = blo->productRange()->numBlocks();
32  int cols = blo->productDomain()->numBlocks();
33 
34  TEUCHOS_ASSERT(rows == cols);
35 
36  // get diagonal blocks
37  std::vector<LinearOp> invDiag;
38  invOpsStrategy_->getInvD(blo, state, invDiag);
39  TEUCHOS_ASSERT(rows == (int)invDiag.size());
40 
41  return createBlockDiagonalInverseOp(blo, invDiag, "Jacobi");
42 }
43 
45 void JacobiPreconditionerFactory::initializeFromParameterList(const Teuchos::ParameterList& pl)
46 #if 0
47 {
48  RCP<const InverseLibrary> invLib = getInverseLibrary();
49 
50  // get string specifying inverse
51  std::string invStr = pl.get<std::string>("Inverse Type");
52 #if defined(Teko_ENABLE_Amesos)
53  if(invStr=="") invStr = "Amesos";
54 #elif defined(Teko_ENABLE_Amesos2)
55  if(invStr=="") invStr= "Amesos2";
56 #endif
57 
58  // based on parameter type build a strategy
59  invOpsStrategy_ = rcp(new InvFactoryDiagStrategy(invLib->getInverseFactory(invStr)));
60 }
61 #endif
62 {
63  Teko_DEBUG_SCOPE("JacobiPreconditionerFactory::initializeFromParameterList", 10);
64  Teko_DEBUG_MSG_BEGIN(9);
65  DEBUG_STREAM << "Parameter list: " << std::endl;
66  pl.print(DEBUG_STREAM);
67  Teko_DEBUG_MSG_END();
68 
69  const std::string inverse_type = "Inverse Type";
70  const std::string preconditioner_type = "Preconditioner Type";
71  std::vector<RCP<InverseFactory> > inverses;
72  std::vector<RCP<InverseFactory> > preconditioners;
73 
74  RCP<const InverseLibrary> invLib = getInverseLibrary();
75 
76  // get string specifying default inverse
77  std::string invStr = "";
78 #if defined(Teko_ENABLE_Amesos)
79  invStr = "Amesos";
80 #elif defined(Teko_ENABLE_Amesos2)
81  invStr = "Amesos2";
82 #endif
83 
84  std::string precStr = "None";
85  if (pl.isParameter(inverse_type)) invStr = pl.get<std::string>(inverse_type);
86  RCP<InverseFactory> defaultPrec;
87  if (precStr != "None") defaultPrec = invLib->getInverseFactory(precStr);
88 
89  Teko_DEBUG_MSG("JacobiPrecFact: Building default inverse \"" << invStr << "\"", 5);
90  RCP<InverseFactory> defaultInverse = invLib->getInverseFactory(invStr);
91 
92  // now check individual solvers
93  Teuchos::ParameterList::ConstIterator itr;
94  for (itr = pl.begin(); itr != pl.end(); ++itr) {
95  std::string fieldName = itr->first;
96  Teko_DEBUG_MSG("JacobiPrecFact: checking fieldName = \"" << fieldName << "\"", 9);
97 
98  // figure out what the integer is
99  if (fieldName.compare(0, inverse_type.length(), inverse_type) == 0 &&
100  fieldName != inverse_type) {
101  int position = -1;
102  std::string inverse, type;
103 
104  // figure out position
105  std::stringstream ss(fieldName);
106  ss >> inverse >> type >> position;
107 
108  if (position <= 0) {
109  Teko_DEBUG_MSG("Jacobi \"Inverse Type\" must be a (strictly) positive integer", 1);
110  }
111 
112  // inserting inverse factory into vector
113  std::string invStr2 = pl.get<std::string>(fieldName);
114  Teko_DEBUG_MSG("JacobiPrecFact: Building inverse " << position << " \"" << invStr2 << "\"",
115  5);
116  if (position > (int)inverses.size()) {
117  inverses.resize(position, defaultInverse);
118  inverses[position - 1] = invLib->getInverseFactory(invStr2);
119  } else
120  inverses[position - 1] = invLib->getInverseFactory(invStr2);
121  } else if (fieldName.compare(0, preconditioner_type.length(), preconditioner_type) == 0 &&
122  fieldName != preconditioner_type) {
123  int position = -1;
124  std::string preconditioner, type;
125 
126  // figure out position
127  std::stringstream ss(fieldName);
128  ss >> preconditioner >> type >> position;
129 
130  if (position <= 0) {
131  Teko_DEBUG_MSG("Jacobi \"Preconditioner Type\" must be a (strictly) positive integer", 1);
132  }
133 
134  // inserting preconditioner factory into vector
135  std::string precStr2 = pl.get<std::string>(fieldName);
136  Teko_DEBUG_MSG(
137  "JacobiPrecFact: Building preconditioner " << position << " \"" << precStr2 << "\"", 5);
138  if (position > (int)preconditioners.size()) {
139  preconditioners.resize(position, defaultPrec);
140  preconditioners[position - 1] = invLib->getInverseFactory(precStr2);
141  } else
142  preconditioners[position - 1] = invLib->getInverseFactory(precStr2);
143  }
144  }
145 
146  // use default inverse
147  if (inverses.size() == 0) inverses.push_back(defaultInverse);
148 
149  // based on parameter type build a strategy
150  invOpsStrategy_ =
151  rcp(new InvFactoryDiagStrategy(inverses, preconditioners, defaultInverse, defaultPrec));
152 }
153 
154 } // namespace Teko
An implementation of a state object for block preconditioners.
Teuchos::RCP< const BlockInvDiagonalStrategy > invOpsStrategy_
some members
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl)
Initialize from a parameter list.
LinearOp buildPreconditionerOperator(BlockedLinearOp &blo, BlockPreconditionerState &state) const
Create the Jacobi preconditioner operator.