Teko Version of the Day
Loading...
Searching...
No Matches
Teko_GaussSeidelPreconditionerFactory.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_GaussSeidelPreconditionerFactory.hpp"
11
12#include "Teko_BlockUpperTriInverseOp.hpp"
13#include "Teko_BlockLowerTriInverseOp.hpp"
14#include "Teuchos_CompilerCodeTweakMacros.hpp"
15
16using Teuchos::rcp;
17using Teuchos::RCP;
18
19namespace Teko {
20
22 const LinearOp& invD0,
23 const LinearOp& invD1)
24 : invOpsStrategy_(rcp(new StaticInvDiagStrategy(invD0, invD1))), solveType_(solveType) {}
25
27 TriSolveType solveType, const RCP<const BlockInvDiagonalStrategy>& strategy)
28 : invOpsStrategy_(strategy), solveType_(solveType) {}
29
32
34 BlockedLinearOp& blo, BlockPreconditionerState& state) const {
35 int rows = blockRowCount(blo);
36 int cols = blockColCount(blo);
37
38 TEUCHOS_ASSERT(rows == cols);
39
40 // get diagonal blocks
41 std::vector<LinearOp> invDiag;
42 invOpsStrategy_->getInvD(blo, state, invDiag);
43 TEUCHOS_ASSERT(rows == (int)invDiag.size());
44
45 if (solveType_ == GS_UseUpperTriangle) {
46 // create a blocked linear operator
47 BlockedLinearOp U = getUpperTriBlocks(blo);
48
49 return createBlockUpperTriInverseOp(U, invDiag, "Gauss Seidel");
50 } else if (solveType_ == GS_UseLowerTriangle) {
51 // create a blocked linear operator
52 BlockedLinearOp L = getLowerTriBlocks(blo);
53
54 return createBlockLowerTriInverseOp(L, invDiag, "Gauss Seidel");
55 }
56
57 TEUCHOS_ASSERT(false); // we should never make it here!
58 TEUCHOS_UNREACHABLE_RETURN(Teuchos::null);
59}
60
63 const Teuchos::ParameterList& pl) {
64 Teko_DEBUG_SCOPE("GaussSeidelPreconditionerFactory::initializeFromParameterList", 10);
65 Teko_DEBUG_MSG_BEGIN(9);
66 DEBUG_STREAM << "Parameter list: " << std::endl;
67 pl.print(DEBUG_STREAM);
68 Teko_DEBUG_MSG_END();
69
70 const std::string inverse_type = "Inverse Type";
71 const std::string preconditioner_type = "Preconditioner Type";
72 std::vector<RCP<InverseFactory> > inverses;
73 std::vector<RCP<InverseFactory> > preconditioners;
74
75 RCP<const InverseLibrary> invLib = getInverseLibrary();
76
77 // get string specifying default inverse
78 std::string invStr = "";
79#if defined(Teko_ENABLE_Amesos)
80 invStr = "Amesos";
81#elif defined(Teko_ENABLE_Amesos2)
82 invStr = "Amesos2";
83#endif
84 std::string precStr = "None";
85 if (pl.isParameter(inverse_type)) invStr = pl.get<std::string>(inverse_type);
86 if (pl.isParameter(preconditioner_type)) precStr = pl.get<std::string>(preconditioner_type);
87 if (pl.isParameter("Use Upper Triangle"))
88 solveType_ = pl.get<bool>("Use Upper Triangle") ? GS_UseUpperTriangle : GS_UseLowerTriangle;
89
90 Teko_DEBUG_MSG("GSPrecFact: Building default inverse \"" << invStr << "\"", 5);
91 RCP<InverseFactory> defaultInverse = invLib->getInverseFactory(invStr);
92 RCP<InverseFactory> defaultPrec;
93 if (precStr != "None") defaultPrec = invLib->getInverseFactory(precStr);
94
95 // now check individual solvers
96 Teuchos::ParameterList::ConstIterator itr;
97 for (itr = pl.begin(); itr != pl.end(); ++itr) {
98 std::string fieldName = itr->first;
99 Teko_DEBUG_MSG("GSPrecFact: checking fieldName = \"" << fieldName << "\"", 9);
100
101 // figure out what the integer is
102 if (fieldName.compare(0, inverse_type.length(), inverse_type) == 0 &&
103 fieldName != inverse_type) {
104 int position = -1;
105 std::string inverse, type;
106
107 // figure out position
108 std::stringstream ss(fieldName);
109 ss >> inverse >> type >> position;
110
111 if (position <= 0) {
112 Teko_DEBUG_MSG("Gauss-Seidel \"Inverse Type\" must be a (strictly) positive integer", 1);
113 }
114
115 // inserting inverse factory into vector
116 std::string invStr2 = pl.get<std::string>(fieldName);
117 Teko_DEBUG_MSG("GSPrecFact: Building inverse " << position << " \"" << invStr2 << "\"", 5);
118 if (position > (int)inverses.size()) {
119 inverses.resize(position, defaultInverse);
120 inverses[position - 1] = invLib->getInverseFactory(invStr2);
121 } else
122 inverses[position - 1] = invLib->getInverseFactory(invStr2);
123 } else if (fieldName.compare(0, preconditioner_type.length(), preconditioner_type) == 0 &&
124 fieldName != preconditioner_type) {
125 int position = -1;
126 std::string preconditioner, type;
127
128 // figure out position
129 std::stringstream ss(fieldName);
130 ss >> preconditioner >> type >> position;
131
132 if (position <= 0) {
133 Teko_DEBUG_MSG("Gauss-Seidel \"Preconditioner Type\" must be a (strictly) positive integer",
134 1);
135 }
136
137 // inserting preconditioner factory into vector
138 std::string precStr2 = pl.get<std::string>(fieldName);
139 Teko_DEBUG_MSG(
140 "GSPrecFact: Building preconditioner " << position << " \"" << precStr2 << "\"", 5);
141 if (position > (int)preconditioners.size()) {
142 preconditioners.resize(position, defaultPrec);
143 preconditioners[position - 1] = invLib->getInverseFactory(precStr2);
144 } else
145 preconditioners[position - 1] = invLib->getInverseFactory(precStr2);
146 }
147 }
148
149 // use default inverse
150 if (inverses.size() == 0) inverses.push_back(defaultInverse);
151
152 // based on parameter type build a strategy
154 rcp(new InvFactoryDiagStrategy(inverses, preconditioners, defaultInverse, defaultPrec));
155}
156
157} // namespace Teko
An implementation of a state object for block preconditioners.
LinearOp buildPreconditionerOperator(BlockedLinearOp &blo, BlockPreconditionerState &state) const
Create the Gauss-Seidel preconditioner operator.
Teuchos::RCP< const BlockInvDiagonalStrategy > invOpsStrategy_
some members
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl)
Initialize from a parameter list.
Teuchos::RCP< const InverseLibrary > getInverseLibrary() const
Get the inverse library used by this preconditioner factory.