10 #include "Teko_JacobiPreconditionerFactory.hpp"
11 #include "Teko_BlockDiagonalInverseOp.hpp"
18 const LinearOp& invD1)
22 const RCP<const BlockInvDiagonalStrategy>& strategy)
23 : invOpsStrategy_(strategy) {}
31 int rows = blo->productRange()->numBlocks();
32 int cols = blo->productDomain()->numBlocks();
34 TEUCHOS_ASSERT(rows == cols);
37 std::vector<LinearOp> invDiag;
39 TEUCHOS_ASSERT(rows == (
int)invDiag.size());
41 return createBlockDiagonalInverseOp(blo, invDiag,
"Jacobi");
48 RCP<const InverseLibrary> invLib = getInverseLibrary();
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";
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);
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;
74 RCP<const InverseLibrary> invLib = getInverseLibrary();
77 std::string invStr =
"";
78 #if defined(Teko_ENABLE_Amesos)
80 #elif defined(Teko_ENABLE_Amesos2)
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);
89 Teko_DEBUG_MSG(
"JacobiPrecFact: Building default inverse \"" << invStr <<
"\"", 5);
90 RCP<InverseFactory> defaultInverse = invLib->getInverseFactory(invStr);
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);
99 if (fieldName.compare(0, inverse_type.length(), inverse_type) == 0 &&
100 fieldName != inverse_type) {
102 std::string inverse, type;
105 std::stringstream ss(fieldName);
106 ss >> inverse >> type >> position;
109 Teko_DEBUG_MSG(
"Jacobi \"Inverse Type\" must be a (strictly) positive integer", 1);
113 std::string invStr2 = pl.get<std::string>(fieldName);
114 Teko_DEBUG_MSG(
"JacobiPrecFact: Building inverse " << position <<
" \"" << invStr2 <<
"\"",
116 if (position > (
int)inverses.size()) {
117 inverses.resize(position, defaultInverse);
118 inverses[position - 1] = invLib->getInverseFactory(invStr2);
120 inverses[position - 1] = invLib->getInverseFactory(invStr2);
121 }
else if (fieldName.compare(0, preconditioner_type.length(), preconditioner_type) == 0 &&
122 fieldName != preconditioner_type) {
124 std::string preconditioner, type;
127 std::stringstream ss(fieldName);
128 ss >> preconditioner >> type >> position;
131 Teko_DEBUG_MSG(
"Jacobi \"Preconditioner Type\" must be a (strictly) positive integer", 1);
135 std::string precStr2 = pl.get<std::string>(fieldName);
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);
142 preconditioners[position - 1] = invLib->getInverseFactory(precStr2);
147 if (inverses.size() == 0) inverses.push_back(defaultInverse);
151 rcp(
new InvFactoryDiagStrategy(inverses, preconditioners, defaultInverse, defaultPrec));
JacobiPreconditionerFactory()
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.