A strategy that takes a single inverse factory and uses that for all inverses. If no mass matrix is passed in the diagonal of the 1,1 block is used. More...
#include <Teko_LSCSIMPLECStrategy.hpp>

Public Member Functions | |
Constructors | |
| LSCSIMPLECStrategy () | |
| virtual void | buildState (BlockedLinearOp &A, BlockPreconditionerState &state) const |
| Functions inherited from LSCStrategy. | |
| virtual LinearOp | getInvBQBt (const BlockedLinearOp &A, BlockPreconditionerState &state) const |
| virtual LinearOp | getInvBHBt (const BlockedLinearOp &A, BlockPreconditionerState &state) const |
| virtual LinearOp | getInvF (const BlockedLinearOp &A, BlockPreconditionerState &state) const |
| virtual LinearOp | getOuterStabilization (const BlockedLinearOp &A, BlockPreconditionerState &state) const |
| virtual LinearOp | getInnerStabilization (const BlockedLinearOp &A, BlockPreconditionerState &state) const |
| virtual LinearOp | getInvMass (const BlockedLinearOp &A, BlockPreconditionerState &state) const |
| virtual LinearOp | getHScaling (const BlockedLinearOp &A, BlockPreconditionerState &state) const |
| virtual bool | useFullLDU () const |
| virtual void | initializeFromParameterList (const Teuchos::ParameterList &pl, const InverseLibrary &invLib) |
| Initialize from a parameter list. | |
| virtual void | initializeState (const BlockedLinearOp &A, LSCPrecondState *state) const |
| Initialize the state object using this blocked linear operator. | |
| void | computeInverses (const BlockedLinearOp &A, LSCPrecondState *state) const |
| virtual void | setUseFullLDU (bool val) |
| Set to true to use the Full LDU decomposition, false otherwise. | |
| virtual void | setSymmetric (bool isSymmetric) |
A strategy that takes a single inverse factory and uses that for all inverses. If no mass matrix is passed in the diagonal of the 1,1 block is used.
A strategy that takes a single inverse factory and uses that for all inverses. Optionally the mass matrix can be passed in, if it is the diagonal is extracted and that is used to form the inverse approximation.
Definition at line 66 of file Teko_LSCSIMPLECStrategy.hpp.
| void Teko::NS::LSCSIMPLECStrategy::buildState | ( | BlockedLinearOp & | A, |
| BlockPreconditionerState & | state | ||
| ) | const [virtual] |
Functions inherited from LSCStrategy.
This informs the strategy object to build the state associated with this operator.
| [in] | A | The linear operator to be preconditioned by LSC. |
| [in] | state | State object for storying reusable information about the operator A. |
Implements Teko::NS::LSCStrategy.
Definition at line 91 of file Teko_LSCSIMPLECStrategy.cpp.
| LinearOp Teko::NS::LSCSIMPLECStrategy::getInvBQBt | ( | const BlockedLinearOp & | A, |
| BlockPreconditionerState & | state | ||
| ) | const [virtual] |
Get the inverse of
.
| [in] | A | The linear operator to be preconditioned by LSC. |
| [in] | state | State object for storying reusable information about the operator A. |
. Implements Teko::NS::LSCStrategy.
Definition at line 127 of file Teko_LSCSIMPLECStrategy.cpp.
| LinearOp Teko::NS::LSCSIMPLECStrategy::getInvBHBt | ( | const BlockedLinearOp & | A, |
| BlockPreconditionerState & | state | ||
| ) | const [virtual] |
Get the inverse of
.
| [in] | A | The linear operator to be preconditioned by LSC. |
| [in] | state | State object for storying reusable information about the operator A. |
. Implements Teko::NS::LSCStrategy.
Definition at line 132 of file Teko_LSCSIMPLECStrategy.cpp.
| LinearOp Teko::NS::LSCSIMPLECStrategy::getInvF | ( | const BlockedLinearOp & | A, |
| BlockPreconditionerState & | state | ||
| ) | const [virtual] |
Get the inverse of the
block.
| [in] | A | The linear operator to be preconditioned by LSC. |
| [in] | state | State object for storying reusable information about the operator A. |
. Implements Teko::NS::LSCStrategy.
Definition at line 137 of file Teko_LSCSIMPLECStrategy.cpp.
| virtual LinearOp Teko::NS::LSCSIMPLECStrategy::getOuterStabilization | ( | const BlockedLinearOp & | A, |
| BlockPreconditionerState & | state | ||
| ) | const [inline, virtual] |
Get the inverse for stabilizing the whole schur complement approximation. Returns
.
| [in] | A | The linear operator to be preconditioned by LSC. |
| [in] | state | State object for storying reusable information about the operator A. |
Implements Teko::NS::LSCStrategy.
Definition at line 126 of file Teko_LSCSIMPLECStrategy.hpp.
| LinearOp Teko::NS::LSCSIMPLECStrategy::getInnerStabilization | ( | const BlockedLinearOp & | A, |
| BlockPreconditionerState & | state | ||
| ) | const [virtual] |
Get the inverse to stablized stabilizing the Schur complement approximation using a placement on the ``inside''. That is what is the value for
. This quantity may be null.
| [in] | A | The linear operator to be preconditioned by LSC. |
| [in] | state | State object for storying reusable information about the operator A. |
Implements Teko::NS::LSCStrategy.
Definition at line 142 of file Teko_LSCSIMPLECStrategy.cpp.
| LinearOp Teko::NS::LSCSIMPLECStrategy::getInvMass | ( | const BlockedLinearOp & | A, |
| BlockPreconditionerState & | state | ||
| ) | const [virtual] |
Get the inverse mass matrix.
| [in] | A | The linear operator to be preconditioned by LSC. |
| [in] | state | State object for storying reusable information about the operator A. |
. Implements Teko::NS::LSCStrategy.
Definition at line 153 of file Teko_LSCSIMPLECStrategy.cpp.
| LinearOp Teko::NS::LSCSIMPLECStrategy::getHScaling | ( | const BlockedLinearOp & | A, |
| BlockPreconditionerState & | state | ||
| ) | const [virtual] |
Get the
scaling matrix.
| [in] | A | The linear operator to be preconditioned by LSC. |
| [in] | state | State object for storying reusable information about the operator A. |
scaling matrix. Implements Teko::NS::LSCStrategy.
Definition at line 162 of file Teko_LSCSIMPLECStrategy.cpp.
| virtual bool Teko::NS::LSCSIMPLECStrategy::useFullLDU | ( | ) | const [inline, virtual] |
Should the approximation of the inverse use a full LDU decomposition, or is a upper triangular approximation sufficient.
Implements Teko::NS::LSCStrategy.
Definition at line 157 of file Teko_LSCSIMPLECStrategy.hpp.
| void Teko::NS::LSCSIMPLECStrategy::initializeFromParameterList | ( | const Teuchos::ParameterList & | pl, |
| const InverseLibrary & | invLib | ||
| ) | [virtual] |
Initialize from a parameter list.
Reimplemented from Teko::NS::LSCStrategy.
Definition at line 239 of file Teko_LSCSIMPLECStrategy.cpp.
| void Teko::NS::LSCSIMPLECStrategy::initializeState | ( | const BlockedLinearOp & | A, |
| LSCPrecondState * | state | ||
| ) | const [virtual] |
Initialize the state object using this blocked linear operator.
Definition at line 168 of file Teko_LSCSIMPLECStrategy.cpp.
| void Teko::NS::LSCSIMPLECStrategy::computeInverses | ( | const BlockedLinearOp & | A, |
| LSCPrecondState * | state | ||
| ) | const |
Compute the inverses required for the LSC Schur complement
Definition at line 199 of file Teko_LSCSIMPLECStrategy.cpp.
| virtual void Teko::NS::LSCSIMPLECStrategy::setUseFullLDU | ( | bool | val | ) | [inline, virtual] |
Set to true to use the Full LDU decomposition, false otherwise.
Definition at line 174 of file Teko_LSCSIMPLECStrategy.hpp.
| virtual void Teko::NS::LSCSIMPLECStrategy::setSymmetric | ( | bool | isSymmetric | ) | [inline, virtual] |
Tell strategy that this operator is supposed to be symmetric. Behavior of LSC is slightly different for non-symmetric case.
| [in] | isSymmetric | Is this operator symmetric? |
Implements Teko::NS::LSCStrategy.
Definition at line 176 of file Teko_LSCSIMPLECStrategy.hpp.
1.7.6.1