Sparse matrix vector multiplication for augmented Lagrangian-based preconditioners. More...
#include <Teko_ALOperator.hpp>
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 MappingStrategy > | getMapStrategy () 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_ |
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:
The equivalent augmented Lagrangian formulation is:
or
Here can be take as the diagonal of the pressure mass matrix and is a positive number.
This class implements the matrix vector product with .
Definition at line 83 of file Teko_ALOperator.hpp.
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.
[in] | vars | Vector of vectors of global ids specifying how the operator is to be blocked. |
[in] | content | Operator to be blocked |
[in] | pressureMassMatrix | Pressure mass matrix |
[in] | gamma | Augmentation parameter |
[in] | label | Label for name the operator |
Definition at line 34 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.
[in] | vars | Vector of vectors of global ids specifying how the operator is to be blocked. |
[in] | content | Operator to be blocked |
[in] | gamma | Augmentation parameter |
[in] | label | Name of the operator |
Definition at line 46 of file Teko_ALOperator.cpp.
void Teko::NS::ALOperator::setPressureMassMatrix | ( | LinearOp | pressureMassMatrix | ) |
Set the pressure mass matrix.
[in] | pressureMassMatrix | Pressure mass matrix. |
Definition at line 58 of file Teko_ALOperator.cpp.
|
inline |
Definition at line 143 of file Teko_ALOperator.hpp.
void Teko::NS::ALOperator::setGamma | ( | double | gamma | ) |
Set gamma.
[in] | gamma | Augmentation parameter. |
Definition at line 67 of file Teko_ALOperator.cpp.
|
inline |
Definition at line 161 of file Teko_ALOperator.hpp.
void Teko::NS::ALOperator::augmentRHS | ( | const Epetra_MultiVector & | b, |
Epetra_MultiVector & | bAugmented | ||
) |
[in] | b | Right-hand side. |
[out] | bAugmented | Augmented right-hand side. |
Definition at line 215 of file Teko_ALOperator.cpp.
|
inline |
Definition at line 179 of file Teko_ALOperator.hpp.
|
inlinevirtual |
Force a rebuild of the blocked operator from the stored content operator.
Reimplemented from Teko::Epetra::BlockedEpetraOperator.
Definition at line 189 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.
[in] | i | Row index. |
[in] | j | Column index. |
Definition at line 74 of file Teko_ALOperator.cpp.
|
protected |
Check dimension. Only implemente for 2D and 3D problems.
Definition at line 84 of file Teko_ALOperator.cpp.
|
protected |
Build AL operator.
Definition at line 91 of file Teko_ALOperator.cpp.
|
protected |
AL operator.
Definition at line 209 of file Teko_ALOperator.hpp.
|
protected |
Operator for augmenting the right-hand side.
Definition at line 214 of file Teko_ALOperator.hpp.
|
protected |
Pressure mass matrix and its inverse.
Definition at line 219 of file Teko_ALOperator.hpp.
|
protected |
Inverse of the pressure mass matrix.
Definition at line 224 of file Teko_ALOperator.hpp.
|
protected |
Augmentation parameter.
Definition at line 229 of file Teko_ALOperator.hpp.
|
protected |
Dimension of the problem.
Definition at line 234 of file Teko_ALOperator.hpp.
|
protected |
Number of block rows.
Definition at line 239 of file Teko_ALOperator.hpp.