Teko  Version of the Day
 All Classes Files Functions Variables Pages
Teko_LSCStrategy.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_LSCStrategy_hpp__
11 #define __Teko_LSCStrategy_hpp__
12 
13 #include "Teuchos_RCP.hpp"
14 
15 #include "Thyra_LinearOpBase.hpp"
16 
17 #include "Teko_Utilities.hpp"
18 #include "Teko_InverseFactory.hpp"
19 #include "Teko_BlockPreconditionerFactory.hpp"
20 
21 namespace Teko {
22 namespace NS {
23 
24 class LSCPrecondState; // forward declaration
25 
82 class LSCStrategy {
83  public:
84  virtual ~LSCStrategy() {}
85 
93  virtual void buildState(BlockedLinearOp& A, BlockPreconditionerState& state) const = 0;
94 
103  virtual LinearOp getInvBQBt(const BlockedLinearOp& A, BlockPreconditionerState& state) const = 0;
104 
113  virtual LinearOp getInvBHBt(const BlockedLinearOp& A, BlockPreconditionerState& state) const = 0;
114 
123  virtual LinearOp getInvF(const BlockedLinearOp& A, BlockPreconditionerState& state) const = 0;
124 
125 #if 0
126 
134  virtual LinearOp getInvAlphaD(const BlockedLinearOp & A,BlockPreconditionerState & state) const = 0;
135 #endif
136 
148  virtual LinearOp getOuterStabilization(const BlockedLinearOp& A,
149  BlockPreconditionerState& state) const = 0;
150 
161  virtual LinearOp getInnerStabilization(const BlockedLinearOp& A,
162  BlockPreconditionerState& state) const = 0;
163 
172  virtual LinearOp getInvMass(const BlockedLinearOp& A, BlockPreconditionerState& state) const = 0;
173 
182  virtual LinearOp getHScaling(const BlockedLinearOp& A, BlockPreconditionerState& state) const = 0;
183 
190  virtual bool useFullLDU() const = 0;
191 
197  virtual void setSymmetric(bool isSymmetric) = 0;
198 
200  virtual void initializeFromParameterList(const Teuchos::ParameterList& /* pl */,
201  const InverseLibrary& /* invLib */) {}
202 
204  virtual Teuchos::RCP<Teuchos::ParameterList> getRequestedParameters() const {
205  return Teuchos::null;
206  }
207 
209  virtual bool updateRequestedParameters(const Teuchos::ParameterList& /* pl */) { return true; }
210 
212  void setRequestHandler(const Teuchos::RCP<RequestHandler>& rh) { requestHandler_ = rh; }
213 
215  Teuchos::RCP<RequestHandler> getRequestHandler() const { return requestHandler_; }
216 
217  private:
218  Teuchos::RCP<RequestHandler> requestHandler_;
219 };
220 
221 } // end namespace NS
222 } // end namespace Teko
223 
224 #endif
virtual LinearOp getOuterStabilization(const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
virtual Teuchos::RCP< Teuchos::ParameterList > getRequestedParameters() const
For assiting in construction of the preconditioner.
virtual LinearOp getInvBQBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
virtual void buildState(BlockedLinearOp &A, BlockPreconditionerState &state) const =0
virtual LinearOp getInvMass(const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
void setRequestHandler(const Teuchos::RCP< RequestHandler > &rh)
This method sets the request handler for this object.
virtual bool useFullLDU() const =0
An implementation of a state object for block preconditioners.
virtual void initializeFromParameterList(const Teuchos::ParameterList &, const InverseLibrary &)
Initialize from a parameter list.
virtual void setSymmetric(bool isSymmetric)=0
virtual LinearOp getInvBHBt(const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
Strategy for driving LSCPreconditionerFactory.
virtual bool updateRequestedParameters(const Teuchos::ParameterList &)
For assiting in construction of the preconditioner.
Teuchos::RCP< RequestHandler > getRequestHandler() const
This method gets the request handler uses by this object.
virtual LinearOp getInnerStabilization(const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
virtual LinearOp getInvF(const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
virtual LinearOp getHScaling(const BlockedLinearOp &A, BlockPreconditionerState &state) const =0