Teko  Version of the Day
 All Classes Files Functions Variables Pages
Public Member Functions | List of all members
Teko::NS::LSCStrategy Class Referenceabstract

Strategy for driving LSCPreconditionerFactory. More...

#include <Teko_LSCStrategy.hpp>

Inheritance diagram for Teko::NS::LSCStrategy:
Inheritance graph
[legend]

Public Member Functions

virtual void buildState (BlockedLinearOp &A, BlockPreconditionerState &state) const =0
 
virtual LinearOp getInvBQBt (const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
 
virtual LinearOp getInvBHBt (const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
 
virtual LinearOp getInvF (const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
 
virtual LinearOp getOuterStabilization (const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
 
virtual LinearOp getInnerStabilization (const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
 
virtual LinearOp getInvMass (const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
 
virtual LinearOp getHScaling (const BlockedLinearOp &A, BlockPreconditionerState &state) const =0
 
virtual bool useFullLDU () const =0
 
virtual void setSymmetric (bool isSymmetric)=0
 
virtual void initializeFromParameterList (const Teuchos::ParameterList &, const InverseLibrary &)
 Initialize from a parameter list. More...
 
virtual Teuchos::RCP
< Teuchos::ParameterList > 
getRequestedParameters () const
 For assiting in construction of the preconditioner. More...
 
virtual bool updateRequestedParameters (const Teuchos::ParameterList &)
 For assiting in construction of the preconditioner. More...
 
void setRequestHandler (const Teuchos::RCP< RequestHandler > &rh)
 This method sets the request handler for this object. More...
 
Teuchos::RCP< RequestHandlergetRequestHandler () const
 This method gets the request handler uses by this object. More...
 

Detailed Description

Strategy for driving LSCPreconditionerFactory.

Strategy for driving the LSCPreconditionerFactory. This class provides all the pieces required by the LSC preconditioner. The intent is that the user can overide them and build there own implementation. Though a fairly substantial implementation is provided in InvLSCStrategy.

The basics of this method can be found in

[1] Elman, Howle, Shadid, Silvester, and Tuminaro, "Least Squares Preconditioners for Stabilized Discretizations of the Navier-Stokes Euqations," SISC-2007.

[2] Elman, and Tuminaro, "Boundary Conditions in Approximate Commutator Preconditioners for the Navier-Stokes Equations," In press (8/2009)?

The Least Squares Commuator preconditioner provides a (nearly) Algebraic approximation of the Schur complement of the (Navier-)Stokes system

$ A = \left[\begin{array}{cc} F & B^T \\ B & C \end{array}\right] $

The approximation to the Schur complement is

$ C - B F^{-1} B^T \approx (B \hat{Q}_u^{-1} B^T - \gamma C)^{-1} (B \hat{Q}_u^{-1} F H B^T+C_I) (B H B^T - \gamma C)^{-1} + C_O $.

Where $\hat{Q}_u$ is typically a diagonal approximation of the mass matrix, and $H$ is an appropriate diagonal scaling matrix (see [2] for details). The scalars $\alpha$ and $\gamma$ are chosen to stabilize an unstable discretization (for the case of $C\neq 0$). If the system is stable then they can be set to $0$ (see [1] for more details).

In order to approximate $A$ two decompositions can be chosen, a full LU decomposition and a purely upper triangular version. A full LU decomposition requires that the velocity convection-diffusion operator ( $F$) is inverted twice, while an upper triangular approximation requires only a single inverse.

The methods of this strategy provide the different pieces. For instance getInvF provides $F^{-1}$. Similarly there are calls to get the inverses of $B \hat{Q}_u^{-1} B^T - \gamma C$, $B \hat{Q}_u^{-1} B^T - \gamma C$, and $\hat{Q}_u^{-1}$ as well as the $H$ operator. All these methods are required by the LSCPreconditionerFactory. Additionally there is a buildState method that is called everytime a preconditiner is (re)constructed. This is to allow for any preprocessing neccessary to be handled.

The final set of methods help construct a LSCStrategy object, they are primarily used by the parameter list construction inteface. They are more advanced and can be ignored by initial implementations of this class.

Definition at line 119 of file Teko_LSCStrategy.hpp.

Member Function Documentation

virtual void Teko::NS::LSCStrategy::buildState ( BlockedLinearOp &  A,
BlockPreconditionerState state 
) const
pure virtual

This informs the strategy object to build the state associated with this operator.

Parameters
[in]AThe linear operator to be preconditioned by LSC.
[in]stateState object for storying reusable information about the operator A.

Implemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

virtual LinearOp Teko::NS::LSCStrategy::getInvBQBt ( const BlockedLinearOp &  A,
BlockPreconditionerState state 
) const
pure virtual

Get the inverse of $B Q_u^{-1} B^T - \gamma C$.

Parameters
[in]AThe linear operator to be preconditioned by LSC.
[in]stateState object for storying reusable information about the operator A.
Returns
An (approximate) inverse of $B Q_u^{-1} B^T - \gamma C$.

Implemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

virtual LinearOp Teko::NS::LSCStrategy::getInvBHBt ( const BlockedLinearOp &  A,
BlockPreconditionerState state 
) const
pure virtual

Get the inverse of $B H B^T - \gamma C$.

Parameters
[in]AThe linear operator to be preconditioned by LSC.
[in]stateState object for storying reusable information about the operator A.
Returns
An (approximate) inverse of $B H B^T - \gamma C$.

Implemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

virtual LinearOp Teko::NS::LSCStrategy::getInvF ( const BlockedLinearOp &  A,
BlockPreconditionerState state 
) const
pure virtual

Get the inverse of the $F$ block.

Parameters
[in]AThe linear operator to be preconditioned by LSC.
[in]stateState object for storying reusable information about the operator A.
Returns
An (approximate) inverse of $F$.

Implemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

virtual LinearOp Teko::NS::LSCStrategy::getOuterStabilization ( const BlockedLinearOp &  A,
BlockPreconditionerState state 
) const
pure virtual

Get the inverse to stablized stabilizing the Schur complement approximation using a placement on the ``outside''. That is what is the value for $C_O$. This quantity may be null.

Parameters
[in]AThe linear operator to be preconditioned by LSC.
[in]stateState object for storying reusable information about the operator A.
Returns
The operator to stabilize the whole Schur complement (originally $\alpha D^{-1} $).

Implemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

virtual LinearOp Teko::NS::LSCStrategy::getInnerStabilization ( const BlockedLinearOp &  A,
BlockPreconditionerState state 
) const
pure virtual

Get the inverse to stablized stabilizing the Schur complement approximation using a placement on the ``inside''. That is what is the value for $C_I$. This quantity may be null.

Parameters
[in]AThe linear operator to be preconditioned by LSC.
[in]stateState object for storying reusable information about the operator A.
Returns
The operator to stabilize the whole Schur complement.

Implemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

virtual LinearOp Teko::NS::LSCStrategy::getInvMass ( const BlockedLinearOp &  A,
BlockPreconditionerState state 
) const
pure virtual

Get the inverse mass matrix.

Parameters
[in]AThe linear operator to be preconditioned by LSC.
[in]stateState object for storying reusable information about the operator A.
Returns
The inverse of the mass matrix $Q_u$.

Implemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

virtual LinearOp Teko::NS::LSCStrategy::getHScaling ( const BlockedLinearOp &  A,
BlockPreconditionerState state 
) const
pure virtual

Get the $H$ scaling matrix.

Parameters
[in]AThe linear operator to be preconditioned by LSC.
[in]stateState object for storying reusable information about the operator A.
Returns
The $H$ scaling matrix.

Implemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

virtual bool Teko::NS::LSCStrategy::useFullLDU ( ) const
pure virtual

Should the approximation of the inverse use a full LDU decomposition, or is a upper triangular approximation sufficient.

Returns
True if the full LDU decomposition should be used, otherwise only an upper triangular version is used.

Implemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

virtual void Teko::NS::LSCStrategy::setSymmetric ( bool  isSymmetric)
pure virtual

Tell strategy that this operator is supposed to be symmetric. Behavior of LSC is slightly different for non-symmetric case.

Parameters
[in]isSymmetricIs this operator symmetric?

Implemented in Teko::NS::LSCSIMPLECStrategy, Teko::NS::InvLSCStrategy, and Teko::NS::PresLaplaceLSCStrategy.

virtual void Teko::NS::LSCStrategy::initializeFromParameterList ( const Teuchos::ParameterList &  ,
const InverseLibrary &   
)
inlinevirtual

Initialize from a parameter list.

Reimplemented in Teko::NS::InvLSCStrategy, Teko::NS::PresLaplaceLSCStrategy, and Teko::NS::LSCSIMPLECStrategy.

Definition at line 237 of file Teko_LSCStrategy.hpp.

virtual Teuchos::RCP<Teuchos::ParameterList> Teko::NS::LSCStrategy::getRequestedParameters ( ) const
inlinevirtual

For assiting in construction of the preconditioner.

Reimplemented in Teko::NS::InvLSCStrategy, and Teko::NS::PresLaplaceLSCStrategy.

Definition at line 241 of file Teko_LSCStrategy.hpp.

virtual bool Teko::NS::LSCStrategy::updateRequestedParameters ( const Teuchos::ParameterList &  )
inlinevirtual

For assiting in construction of the preconditioner.

Reimplemented in Teko::NS::InvLSCStrategy, and Teko::NS::PresLaplaceLSCStrategy.

Definition at line 246 of file Teko_LSCStrategy.hpp.

void Teko::NS::LSCStrategy::setRequestHandler ( const Teuchos::RCP< RequestHandler > &  rh)
inline

This method sets the request handler for this object.

Definition at line 249 of file Teko_LSCStrategy.hpp.

Teuchos::RCP<RequestHandler> Teko::NS::LSCStrategy::getRequestHandler ( ) const
inline

This method gets the request handler uses by this object.

Definition at line 252 of file Teko_LSCStrategy.hpp.


The documentation for this class was generated from the following file: