Teko  Version of the Day
 All Classes Files Functions Variables Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
Teko::NS::ALOperator Class Reference

Sparse matrix vector multiplication for augmented Lagrangian-based preconditioners. More...

#include <Teko_ALOperator.hpp>

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

Public Member Functions

 ALOperator (const std::vector< std::vector< int > > &vars, const Teuchos::RCP< Epetra_Operator > &content, LinearOp pressureMassMatrix, double gamma=0.05, const std::string &label="<ANYM>")
 
 ALOperator (const std::vector< std::vector< int > > &vars, const Teuchos::RCP< Epetra_Operator > &content, double gamma=0.05, const std::string &label="<ANYM>")
 
void setPressureMassMatrix (LinearOp pressureMassMatrix)
 
const LinearOp & getPressureMassMatrix () const
 
void setGamma (double gamma)
 
const double & getGamma () const
 
void augmentRHS (const Epetra_MultiVector &b, Epetra_MultiVector &bAugmented)
 
int getNumberOfBlockRows () const
 
virtual void RebuildOps ()
 
const Teuchos::RCP< const
Epetra_Operator > 
GetBlock (int i, int j) const
 
- Public Member Functions inherited from Teko::Epetra::BlockedEpetraOperator
 BlockedEpetraOperator (const std::vector< std::vector< int > > &vars, const Teuchos::RCP< const Epetra_Operator > &content, const std::string &label="<ANYM>")
 
virtual void SetContent (const std::vector< std::vector< int > > &vars, const Teuchos::RCP< const Epetra_Operator > &content)
 
void Reorder (const BlockReorderManager &brm)
 
void RemoveReording ()
 Remove any reordering on this object. More...
 
virtual void WriteBlocks (const std::string &prefix) const
 
bool testAgainstFullOperator (int count, double tol) const
 Helps perform sanity checks. More...
 
- Public Member Functions inherited from Teko::Epetra::EpetraOperatorWrapper
const RCP< const
Thyra::LinearOpBase< double > > 
getThyraOp () const
 Return the thyra operator associated with this wrapper. More...
 
const RCP< const MappingStrategygetMapStrategy () const
 Get the mapping strategy for this wrapper (translate between Thyra and Epetra) More...
 
virtual int GetBlockRowCount ()
 Get the number of block rows in this operator. More...
 
virtual int GetBlockColCount ()
 Get the number of block columns in this operator. More...
 
Teuchos::RCP< const
Epetra_Operator > 
GetBlock (int i, int j) const
 Grab the i,j block. More...
 

Protected Member Functions

void checkDim (const std::vector< std::vector< int > > &vars)
 
void BuildALOperator ()
 

Protected Attributes

Teuchos::RCP
< Thyra::LinearOpBase< double > > 
alOperator_
 
Teuchos::RCP
< Thyra::LinearOpBase< double > > 
alOperatorRhs_
 
LinearOp pressureMassMatrix_
 
LinearOp invPressureMassMatrix_
 
double gamma_
 
int dim_
 
int numBlockRows_
 

Detailed Description

Sparse matrix vector multiplication for augmented Lagrangian-based preconditioners.

This class implements sparse matrix vector multiplication for augmented Lagrangian-based preconditioners. Details can be found in the following papers:

[1] M. Benzi and M. A. Olshanskii, An Augmented Lagrangian-Based Approach to the Oseen Problem, SIAM J. Scientific Computing, 28 (2006), pp. 2095-2113.

[2] Benzi, M. A. Olshanskii and Z. Wang, Modified Augmented Lagrangian Preconditioners for the Incompressible Navier-Stokes Equations, International Journal for Numerical Methods in Fluids, 66 (2011), pp. 486-508.

Suppose we are solving the following linear system:

$ \left[ \begin{array}{cc} A & B^T \\ B & -C \end{array} \right] \left[ \begin{array}{c} u \\ p \end{array} \right] = \left[ \begin{array}{c} f \\ g \end{array} \right]. $

The equivalent augmented Lagrangian formulation is:

$ \left[ \begin{array}{cc} A + \gamma B^T W^{-1} B & B^T - \gamma B^T W^{-1} C \\ B & -C \end{array} \right] \left[ \begin{array}{c} u \\ p \end{array} \right] = \left[ \begin{array}{c} f + \gamma B^T W^{-1} g \\ g \end{array} \right] $

or

$ \widehat{\mathcal{A}} x = \hat{b}. $

Here $ W $ can be take as the diagonal of the pressure mass matrix and $ \gamma $ is a positive number.

This class implements the matrix vector product with $ \widehat{\mathcal{A}} $.

Definition at line 81 of file Teko_ALOperator.hpp.

Constructor & Destructor Documentation

Teko::NS::ALOperator::ALOperator ( const std::vector< std::vector< int > > &  vars,
const Teuchos::RCP< Epetra_Operator > &  content,
LinearOp  pressureMassMatrix,
double  gamma = 0.05,
const std::string &  label = "<ANYM>" 
)

Build an augmented Lagrangian operator based on a vector of vector of global IDs.

Parameters
[in]varsVector of vectors of global ids specifying how the operator is to be blocked.
[in]contentOperator to be blocked
[in]pressureMassMatrixPressure mass matrix
[in]gammaAugmentation parameter
[in]labelLabel for name the operator

Definition at line 32 of file Teko_ALOperator.cpp.

Teko::NS::ALOperator::ALOperator ( const std::vector< std::vector< int > > &  vars,
const Teuchos::RCP< Epetra_Operator > &  content,
double  gamma = 0.05,
const std::string &  label = "<ANYM>" 
)

Build a modified augmented Lagrangian operator based on a vector of vector of global IDs.

Parameters
[in]varsVector of vectors of global ids specifying how the operator is to be blocked.
[in]contentOperator to be blocked
[in]gammaAugmentation parameter
[in]labelName of the operator

Definition at line 43 of file Teko_ALOperator.cpp.

Member Function Documentation

void Teko::NS::ALOperator::setPressureMassMatrix ( LinearOp  pressureMassMatrix)

Set the pressure mass matrix.

Parameters
[in]pressureMassMatrixPressure mass matrix.

Definition at line 54 of file Teko_ALOperator.cpp.

const LinearOp& Teko::NS::ALOperator::getPressureMassMatrix ( ) const
inline
Returns
Pressure mass matrix that can be used to construct preconditioner.

Definition at line 133 of file Teko_ALOperator.hpp.

void Teko::NS::ALOperator::setGamma ( double  gamma)

Set gamma.

Parameters
[in]gammaAugmentation parameter.

Definition at line 60 of file Teko_ALOperator.cpp.

const double& Teko::NS::ALOperator::getGamma ( ) const
inline
Returns
Gamma Augmentation parameter.

Definition at line 146 of file Teko_ALOperator.hpp.

void Teko::NS::ALOperator::augmentRHS ( const Epetra_MultiVector &  b,
Epetra_MultiVector &  bAugmented 
)
Parameters
[in]bRight-hand side.
[out]bAugmentedAugmented right-hand side.

Definition at line 180 of file Teko_ALOperator.cpp.

int Teko::NS::ALOperator::getNumberOfBlockRows ( ) const
inline
Returns
Number of rows.

Definition at line 159 of file Teko_ALOperator.hpp.

virtual void Teko::NS::ALOperator::RebuildOps ( )
inlinevirtual

Force a rebuild of the blocked operator from the stored content operator.

Reimplemented from Teko::Epetra::BlockedEpetraOperator.

Definition at line 165 of file Teko_ALOperator.hpp.

const Teuchos::RCP< const Epetra_Operator > Teko::NS::ALOperator::GetBlock ( int  i,
int  j 
) const

Get the (i,j) block of the original (non-augmented) operator.

Parameters
[in]iRow index.
[in]jColumn index.

Definition at line 65 of file Teko_ALOperator.cpp.

void Teko::NS::ALOperator::checkDim ( const std::vector< std::vector< int > > &  vars)
protected

Check dimension. Only implemente for 2D and 3D problems.

Definition at line 72 of file Teko_ALOperator.cpp.

void Teko::NS::ALOperator::BuildALOperator ( )
protected

Build AL operator.

Definition at line 77 of file Teko_ALOperator.cpp.

Member Data Documentation

Teuchos::RCP<Thyra::LinearOpBase<double> > Teko::NS::ALOperator::alOperator_
protected

AL operator.

Definition at line 180 of file Teko_ALOperator.hpp.

Teuchos::RCP<Thyra::LinearOpBase<double> > Teko::NS::ALOperator::alOperatorRhs_
protected

Operator for augmenting the right-hand side.

Definition at line 185 of file Teko_ALOperator.hpp.

LinearOp Teko::NS::ALOperator::pressureMassMatrix_
protected

Pressure mass matrix and its inverse.

Definition at line 190 of file Teko_ALOperator.hpp.

LinearOp Teko::NS::ALOperator::invPressureMassMatrix_
protected

Inverse of the pressure mass matrix.

Definition at line 195 of file Teko_ALOperator.hpp.

double Teko::NS::ALOperator::gamma_
protected

Augmentation parameter.

Definition at line 200 of file Teko_ALOperator.hpp.

int Teko::NS::ALOperator::dim_
protected

Dimension of the problem.

Definition at line 205 of file Teko_ALOperator.hpp.

int Teko::NS::ALOperator::numBlockRows_
protected

Number of block rows.

Definition at line 210 of file Teko_ALOperator.hpp.


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