Teko Version of the Day
Loading...
Searching...
No Matches
Teko_InvLSCStrategy.hpp
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#ifndef __Teko_InvLSCStrategy_hpp__
11#define __Teko_InvLSCStrategy_hpp__
12
13#include "Teko_LSCStrategy.hpp"
14
15namespace Teko {
16namespace NS {
17
18class LSCPrecondState; // forward declration
19
29class InvLSCStrategy : public LSCStrategy {
30 public:
32
33 InvLSCStrategy();
34 InvLSCStrategy(const Teuchos::RCP<InverseFactory> &factory, bool rzn = false);
35 InvLSCStrategy(const Teuchos::RCP<InverseFactory> &factory, LinearOp &mass, bool rzn = false);
36 InvLSCStrategy(const Teuchos::RCP<InverseFactory> &invFactF,
37 const Teuchos::RCP<InverseFactory> &invFactS, bool rzn = false);
38 InvLSCStrategy(const Teuchos::RCP<InverseFactory> &invFactF,
39 const Teuchos::RCP<InverseFactory> &invFactS, LinearOp &mass, bool rzn = false);
41
42 virtual ~InvLSCStrategy() {}
43
45
46
54 virtual void buildState(BlockedLinearOp &A, BlockPreconditionerState &state) const;
55
64 virtual LinearOp getInvBQBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
65
74 virtual LinearOp getInvBHBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
75
84 virtual LinearOp getInvF(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
85
94 // virtual LinearOp getInvAlphaD(const BlockedLinearOp & A,BlockPreconditionerState & state)
95 // const;
96 virtual LinearOp getOuterStabilization(const BlockedLinearOp &A,
97 BlockPreconditionerState &state) const;
98 virtual LinearOp getInnerStabilization(const BlockedLinearOp & /* A */,
99 BlockPreconditionerState & /* state */) const {
100 return Teuchos::null;
101 }
102
111 virtual LinearOp getInvMass(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
112
121 virtual LinearOp getHScaling(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
122
129 virtual bool useFullLDU() const { return useFullLDU_; }
130
136 virtual void setSymmetric(bool isSymmetric) { isSymmetric_ = isSymmetric; }
137
139 virtual void initializeFromParameterList(const Teuchos::ParameterList &pl,
140 const InverseLibrary &invLib);
141
143 virtual Teuchos::RCP<Teuchos::ParameterList> getRequestedParameters() const;
144
146 virtual bool updateRequestedParameters(const Teuchos::ParameterList &pl);
148
151 virtual void setPressureStabMatrix(const Teko::LinearOp &psm) { userPresStabMat_ = psm; }
152
154 virtual void initializeState(const BlockedLinearOp &A, LSCPrecondState *state) const;
155
161 void computeInverses(const BlockedLinearOp &A, LSCPrecondState *state) const;
162
163 // //! Initialize the state object using this blocked linear operator
164 // virtual void reinitializeState(const BlockedLinearOp & A,LSCPrecondState * state) const;
165
167 virtual void setEigSolveParam(int sz) { eigSolveParam_ = sz; }
168
170 virtual int getEigSolveParam() { return eigSolveParam_; }
171
173 virtual void setUseFullLDU(bool val) { useFullLDU_ = val; }
174
176 virtual void setRowZeroing(bool val) { rowZeroingNeeded_ = val; }
177
179 virtual void setMassMatrix(const LinearOp &mass) { massMatrix_ = mass; }
180
184 virtual void setHScaling(const LinearOp &hScaling) { hScaling_ = hScaling; }
185
189 virtual void setHScaling(const MultiVector &hScaling) {
190 hScaling_ = buildDiagonal(hScaling, "H");
191 }
192
196 virtual void setWScaling(const MultiVector &wScaling) { wScaling_ = wScaling; }
197
198 protected:
199 LinearOp massMatrix_;
200
201 // how to invert the matrices
202 Teuchos::RCP<InverseFactory> invFactoryF_;
203 Teuchos::RCP<InverseFactory> invFactoryS_;
204
205 // number of power iterations when computing spectral radius
206 int eigSolveParam_;
207
208 // flags for handling various options
209 bool rowZeroingNeeded_;
210 bool useFullLDU_;
211 bool useMass_;
212 bool useLumping_;
213 bool useWScaling_;
214 DiagonalType scaleType_;
215 bool isSymmetric_;
216 bool assumeStable_;
217
218 // operators requested, to be filled by user
219 LinearOp userPresStabMat_;
220 mutable LinearOp hScaling_;
221 MultiVector wScaling_;
222
223 private:
225};
226
227} // end namespace NS
228} // end namespace Teko
229
230#endif
An implementation of a state object for block preconditioners.
A strategy that takes a single inverse factory and uses that for all inverses. If no mass matrix is p...
virtual bool useFullLDU() const
virtual void setHScaling(const MultiVector &hScaling)
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl, const InverseLibrary &invLib)
Initialize from a parameter list.
virtual LinearOp getInvMass(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setSymmetric(bool isSymmetric)
virtual void setMassMatrix(const LinearOp &mass)
set the mass matrix to use in computing the scaling
virtual LinearOp getHScaling(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setUseFullLDU(bool val)
Set to true to use the Full LDU decomposition, false otherwise.
virtual LinearOp getInvF(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setWScaling(const MultiVector &wScaling)
virtual void setEigSolveParam(int sz)
Set the number of power series iterations to use when finding the spectral radius.
virtual void buildState(BlockedLinearOp &A, BlockPreconditionerState &state) const
Functions inherited from LSCStrategy.
virtual LinearOp getOuterStabilization(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual int getEigSolveParam()
Return the number of power series iterations to use when finding the spectral radius.
virtual bool updateRequestedParameters(const Teuchos::ParameterList &pl)
For assiting in construction of the preconditioner.
virtual void initializeState(const BlockedLinearOp &A, LSCPrecondState *state) const
Initialize the state object using this blocked linear operator.
virtual LinearOp getInvBHBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual Teuchos::RCP< Teuchos::ParameterList > getRequestedParameters() const
For assiting in construction of the preconditioner.
void computeInverses(const BlockedLinearOp &A, LSCPrecondState *state) const
virtual void setRowZeroing(bool val)
Set to true to zero the rows of F when computing the spectral radius.
virtual LinearOp getInvBQBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setHScaling(const LinearOp &hScaling)
virtual void setPressureStabMatrix(const Teko::LinearOp &psm)
Preconditioner state for the LSC factory.
Strategy for driving LSCPreconditionerFactory.