10 #ifndef __Teko_InvLSCStrategy_hpp__
11 #define __Teko_InvLSCStrategy_hpp__
13 #include "Teko_LSCStrategy.hpp"
18 class LSCPrecondState;
34 InvLSCStrategy(
const Teuchos::RCP<InverseFactory> &factory,
bool rzn =
false);
35 InvLSCStrategy(
const Teuchos::RCP<InverseFactory> &factory, LinearOp &mass,
bool rzn =
false);
37 const Teuchos::RCP<InverseFactory> &invFactS,
bool rzn =
false);
39 const Teuchos::RCP<InverseFactory> &invFactS, LinearOp &mass,
bool rzn =
false);
100 return Teuchos::null;
136 virtual void setSymmetric(
bool isSymmetric) { isSymmetric_ = isSymmetric; }
140 const InverseLibrary &invLib);
184 virtual void setHScaling(
const LinearOp &hScaling) { hScaling_ = hScaling; }
190 hScaling_ = buildDiagonal(hScaling,
"H");
196 virtual void setWScaling(
const MultiVector &wScaling) { wScaling_ = wScaling; }
199 LinearOp massMatrix_;
202 Teuchos::RCP<InverseFactory> invFactoryF_;
203 Teuchos::RCP<InverseFactory> invFactoryS_;
209 bool rowZeroingNeeded_;
214 DiagonalType scaleType_;
219 LinearOp userPresStabMat_;
220 mutable LinearOp hScaling_;
221 MultiVector wScaling_;
virtual LinearOp getInvF(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setUseFullLDU(bool val)
Set to true to use the Full LDU decomposition, false otherwise.
virtual LinearOp getInnerStabilization(const BlockedLinearOp &, BlockPreconditionerState &) const
virtual void setRowZeroing(bool val)
Set to true to zero the rows of F when computing the spectral radius.
virtual Teuchos::RCP< Teuchos::ParameterList > getRequestedParameters() const
For assiting in construction of the preconditioner.
virtual LinearOp getOuterStabilization(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual void setHScaling(const LinearOp &hScaling)
virtual LinearOp getHScaling(const BlockedLinearOp &A, BlockPreconditionerState &state) const
virtual LinearOp getInvBHBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const
An implementation of a state object for block preconditioners.
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 void initializeState(const BlockedLinearOp &A, LSCPrecondState *state) const
Initialize the state object using this blocked linear operator.
A strategy that takes a single inverse factory and uses that for all inverses. If no mass matrix is p...
virtual LinearOp getInvMass(const BlockedLinearOp &A, BlockPreconditionerState &state) const
void computeInverses(const BlockedLinearOp &A, LSCPrecondState *state) const
Preconditioner state for the LSC factory.
virtual bool updateRequestedParameters(const Teuchos::ParameterList &pl)
For assiting in construction of the preconditioner.
Strategy for driving LSCPreconditionerFactory.
virtual void initializeFromParameterList(const Teuchos::ParameterList &pl, const InverseLibrary &invLib)
Initialize from a parameter list.
virtual void setMassMatrix(const LinearOp &mass)
set the mass matrix to use in computing the scaling
virtual void setHScaling(const MultiVector &hScaling)
virtual int getEigSolveParam()
Return the number of power series iterations to use when finding the spectral radius.
virtual bool useFullLDU() const
virtual void setPressureStabMatrix(const Teko::LinearOp &psm)
virtual void setSymmetric(bool isSymmetric)
virtual void setWScaling(const MultiVector &wScaling)
virtual LinearOp getInvBQBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const