Teko Version of the Day
Loading...
Searching...
No Matches
Teko_PresLaplaceLSCStrategy.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_PresLaplaceLSCStrategy_hpp__
11#define __Teko_PresLaplaceLSCStrategy_hpp__
12
13#include "Teko_LSCStrategy.hpp"
14
15namespace Teko {
16namespace NS {
17
18class LSCPrecondState; // forward declration
19
29class PresLaplaceLSCStrategy : public LSCStrategy {
30 public:
32
33 PresLaplaceLSCStrategy();
34 PresLaplaceLSCStrategy(const Teuchos::RCP<InverseFactory> &factory);
35 PresLaplaceLSCStrategy(const Teuchos::RCP<InverseFactory> &invFactF,
36 const Teuchos::RCP<InverseFactory> &invFactS);
38
39 virtual ~PresLaplaceLSCStrategy() {}
40
42
43
51 virtual void buildState(BlockedLinearOp &A, BlockPreconditionerState &state) const;
52
61 virtual LinearOp getInvBQBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
62
71 virtual LinearOp getInvBHBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
72
81 virtual LinearOp getInvF(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
82
91 // virtual LinearOp getInvAlphaD(const BlockedLinearOp & A,BlockPreconditionerState & state)
92 // const;
93 virtual LinearOp getOuterStabilization(const BlockedLinearOp &A,
94 BlockPreconditionerState &state) const;
95
96 virtual LinearOp getInnerStabilization(const BlockedLinearOp & /* A */,
97 BlockPreconditionerState & /* state */) const {
98 return Teuchos::null;
99 }
100
109 virtual LinearOp getInvMass(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
110
119 virtual LinearOp getHScaling(const BlockedLinearOp &A, BlockPreconditionerState &state) const;
120
127 virtual bool useFullLDU() const { return useFullLDU_; }
128
134 virtual void setSymmetric(bool isSymmetric) { isSymmetric_ = isSymmetric; }
135
137 virtual void initializeFromParameterList(const Teuchos::ParameterList &pl,
138 const InverseLibrary &invLib);
139
141 virtual Teuchos::RCP<Teuchos::ParameterList> getRequestedParameters() const;
142
144 virtual bool updateRequestedParameters(const Teuchos::ParameterList &pl);
146
148 virtual void initializeState(const BlockedLinearOp &A, LSCPrecondState *state) const;
149
155 void computeInverses(const BlockedLinearOp &A, LSCPrecondState *state) const;
156
158 virtual void setEigSolveParam(int sz) { eigSolveParam_ = sz; }
159
161 virtual int getEigSolveParam() { return eigSolveParam_; }
162
164 virtual void setUseFullLDU(bool val) { useFullLDU_ = val; }
165
166 protected:
167 // how to invert the matrices
168 Teuchos::RCP<InverseFactory> invFactoryV_;
169 Teuchos::RCP<InverseFactory> invFactoryP_;
170
171 // flags for handling various options
172 bool isSymmetric_;
173 int eigSolveParam_;
174 bool useFullLDU_;
175
176 // scaling operator parameters
177 bool useMass_;
178 DiagonalType scaleType_;
179
180 private:
182
183 public:
184 // some static functions for determining strings
185
186 static std::string getPressureLaplaceString() { return "Pressure Laplace Operator"; }
187 static std::string getVelocityMassString() { return "Velocity Mass Operator"; }
188};
189
190} // end namespace NS
191} // end namespace Teko
192
193#endif
An implementation of a state object for block preconditioners.
Preconditioner state for the LSC factory.
Strategy for driving LSCPreconditionerFactory.
A strategy that takes a single inverse factory and uses that for all inverses. If no mass matrix is p...
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl, const InverseLibrary &invLib)
Initialize from a parameter list.
virtual int getEigSolveParam()
Return the number of power series iterations to use when finding the spectral radius.
virtual void setSymmetric(bool isSymmetric)
virtual void initializeState(const BlockedLinearOp &A, LSCPrecondState *state) const
Initialize the state object using this blocked linear operator.
virtual LinearOp getHScaling(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual LinearOp getInvMass(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual bool updateRequestedParameters(const Teuchos::ParameterList &pl)
For assiting in construction of the preconditioner.
virtual LinearOp getInvBQBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual Teuchos::RCP< Teuchos::ParameterList > getRequestedParameters() const
For assiting in construction of the preconditioner.
virtual LinearOp getInvF(const BlockedLinearOp &A, BlockPreconditionerState &state) const
void computeInverses(const BlockedLinearOp &A, LSCPrecondState *state) const
virtual LinearOp getOuterStabilization(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void buildState(BlockedLinearOp &A, BlockPreconditionerState &state) const
Functions inherited from LSCStrategy.
virtual void setUseFullLDU(bool val)
Set to true to use the Full LDU decomposition, false otherwise.
virtual void setEigSolveParam(int sz)
Set the number of power series iterations to use when finding the spectral radius.
virtual LinearOp getInvBHBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const