NOX  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Groups Pages
Public Member Functions | Protected Member Functions | Protected Attributes | List of all members
LOCA::BorderedSolver::Bordering Class Reference

Bordered system solver strategy based on bordering. More...

#include <LOCA_BorderedSolver_Bordering.H>

Inheritance diagram for LOCA::BorderedSolver::Bordering:
Inheritance graph
[legend]
Collaboration diagram for LOCA::BorderedSolver::Bordering:
Collaboration graph
[legend]

Public Member Functions

 Bordering (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &solverParams)
 Constructor. More...
 
virtual ~Bordering ()
 Destructor.
 
virtual void setMatrixBlocks (const Teuchos::RCP< const LOCA::BorderedSolver::AbstractOperator > &op, const Teuchos::RCP< const NOX::Abstract::MultiVector > &blockA, const Teuchos::RCP< const LOCA::MultiContinuation::ConstraintInterface > &blockB, const Teuchos::RCP< const NOX::Abstract::MultiVector::DenseMatrix > &blockC)
 Set blocks. More...
 
virtual
NOX::Abstract::Group::ReturnType 
initForSolve ()
 Intialize solver for a solve. More...
 
virtual
NOX::Abstract::Group::ReturnType 
initForTransposeSolve ()
 Intialize solver for a transpose solve. More...
 
virtual
NOX::Abstract::Group::ReturnType 
apply (const NOX::Abstract::MultiVector &X, const NOX::Abstract::MultiVector::DenseMatrix &Y, NOX::Abstract::MultiVector &U, NOX::Abstract::MultiVector::DenseMatrix &V) const
 Computed extended matrix-multivector product. More...
 
virtual
NOX::Abstract::Group::ReturnType 
applyTranspose (const NOX::Abstract::MultiVector &X, const NOX::Abstract::MultiVector::DenseMatrix &Y, NOX::Abstract::MultiVector &U, NOX::Abstract::MultiVector::DenseMatrix &V) const
 Computed extended matrix transpose-multivector product. More...
 
virtual
NOX::Abstract::Group::ReturnType 
applyInverse (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector *F, const NOX::Abstract::MultiVector::DenseMatrix *G, NOX::Abstract::MultiVector &X, NOX::Abstract::MultiVector::DenseMatrix &Y) const
 Solves the extended system as defined above using bordering. More...
 
virtual
NOX::Abstract::Group::ReturnType 
applyInverseTranspose (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector *F, const NOX::Abstract::MultiVector::DenseMatrix *G, NOX::Abstract::MultiVector &X, NOX::Abstract::MultiVector::DenseMatrix &Y) const
 Solves the transpose of the extended system as defined above using bordering. More...
 
- Public Member Functions inherited from LOCA::BorderedSolver::AbstractStrategy
 AbstractStrategy ()
 Constructor.
 
virtual ~AbstractStrategy ()
 Destructor.
 
virtual void setMatrixBlocksMultiVecConstraint (const Teuchos::RCP< const LOCA::BorderedSolver::AbstractOperator > &op, const Teuchos::RCP< const NOX::Abstract::MultiVector > &blockA, const Teuchos::RCP< const NOX::Abstract::MultiVector > &blockB, const Teuchos::RCP< const NOX::Abstract::MultiVector::DenseMatrix > &blockC)
 Set blocks with multivector constraint. More...
 

Protected Member Functions

NOX::Abstract::Group::ReturnType solveFZero (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector *AA, const LOCA::MultiContinuation::ConstraintInterface *BB, const NOX::Abstract::MultiVector::DenseMatrix *CC, const NOX::Abstract::MultiVector::DenseMatrix *G, NOX::Abstract::MultiVector &X, NOX::Abstract::MultiVector::DenseMatrix &Y) const
 Solves the extended system when F = 0.
 
NOX::Abstract::Group::ReturnType solveContiguous (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector *AA, const LOCA::MultiContinuation::ConstraintInterface *BB, const NOX::Abstract::MultiVector::DenseMatrix *CC, std::vector< int > &indexF, std::vector< int > &indexA, const NOX::Abstract::MultiVector *F, const NOX::Abstract::MultiVector::DenseMatrix *G, NOX::Abstract::MultiVector &X, NOX::Abstract::MultiVector::DenseMatrix &Y) const
 Solves the extended system when F and A are contiguous.
 
NOX::Abstract::Group::ReturnType solveFZeroTrans (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector *AA, const NOX::Abstract::MultiVector *BB, const NOX::Abstract::MultiVector::DenseMatrix *CC, const NOX::Abstract::MultiVector::DenseMatrix *G, NOX::Abstract::MultiVector &X, NOX::Abstract::MultiVector::DenseMatrix &Y) const
 Solves the transpose of the extended system when F = 0.
 
NOX::Abstract::Group::ReturnType solveContiguousTrans (Teuchos::ParameterList &params, const NOX::Abstract::MultiVector *AA, const NOX::Abstract::MultiVector *BB, const NOX::Abstract::MultiVector::DenseMatrix *CC, std::vector< int > &indexF, std::vector< int > &indexA, const NOX::Abstract::MultiVector *F, const NOX::Abstract::MultiVector::DenseMatrix *G, NOX::Abstract::MultiVector &X, NOX::Abstract::MultiVector::DenseMatrix &Y) const
 Solves the transpose of the extended system when F and B are contiguous.
 

Protected Attributes

Teuchos::RCP< LOCA::GlobalDataglobalData
 Global data object.
 
Teuchos::RCP
< Teuchos::ParameterList
solverParams
 Solver parameters.
 
Teuchos::RCP< const
LOCA::BorderedSolver::AbstractOperator
op
 Pointer to operator.
 
Teuchos::RCP< const
NOX::Abstract::MultiVector
A
 Pointer to A block.
 
Teuchos::RCP< const
LOCA::MultiContinuation::ConstraintInterface
B
 Pointer to B block.
 
Teuchos::RCP< const
NOX::Abstract::MultiVector::DenseMatrix
C
 Pointer to C block.
 
bool isZeroA
 flag indicating whether A block is zero
 
bool isZeroB
 flag indicating whether B block is zero
 
bool isZeroC
 flag indicating whether C block is zero
 
bool isZeroF
 flag indicating whether F block is zero
 
bool isZeroG
 flag indicating whether G block is zero
 

Detailed Description

Bordered system solver strategy based on bordering.

This class solves the extended system of equations

\[ \begin{bmatrix} J & A \\ B^T & C \end{bmatrix} \begin{bmatrix} X \\ Y \end{bmatrix} = \begin{bmatrix} F \\ G \end{bmatrix} \]

via bordering (block elimination):

\[ \begin{aligned} X_1 &= J^{-1} F \\ X_2 &= J^{-1} A \\ Y &= (C-B^T X_2)^{-1}(G-B^T X_1) \\ X &= X_1 - X_2 Y \end{aligned} \]

It takes advantage of any of the matrix blocks being zero and concatenates $F$ and $A$ into a contiguous multivector to compute $X_1$ and $X_2$ in one block solve.

To solve the transpose of the system, a similar bordering algorithm is implemented. Note however that for the transpose, the constraint object representing $B$ must implement the LOCA::MultiContinuation::ConstraintInterfaceMVDX since $B$ appears on the right-hand-side of a linear system.

Constructor & Destructor Documentation

LOCA::BorderedSolver::Bordering::Bordering ( const Teuchos::RCP< LOCA::GlobalData > &  global_data,
const Teuchos::RCP< LOCA::Parameter::SublistParser > &  topParams,
const Teuchos::RCP< Teuchos::ParameterList > &  solverParams 
)

Constructor.

Parameters
global_data[in] Global data object
topParams[in] Parsed top-level parameter list
solverParams[in] Bordered solver parameters. Currently none are referenced.

Member Function Documentation

NOX::Abstract::Group::ReturnType LOCA::BorderedSolver::Bordering::applyInverse ( Teuchos::ParameterList params,
const NOX::Abstract::MultiVector F,
const NOX::Abstract::MultiVector::DenseMatrix G,
NOX::Abstract::MultiVector X,
NOX::Abstract::MultiVector::DenseMatrix Y 
) const
virtual
NOX::Abstract::Group::ReturnType LOCA::BorderedSolver::Bordering::applyInverseTranspose ( Teuchos::ParameterList params,
const NOX::Abstract::MultiVector F,
const NOX::Abstract::MultiVector::DenseMatrix G,
NOX::Abstract::MultiVector X,
NOX::Abstract::MultiVector::DenseMatrix Y 
) const
virtual
NOX::Abstract::Group::ReturnType LOCA::BorderedSolver::Bordering::applyTranspose ( const NOX::Abstract::MultiVector X,
const NOX::Abstract::MultiVector::DenseMatrix Y,
NOX::Abstract::MultiVector U,
NOX::Abstract::MultiVector::DenseMatrix V 
) const
virtual
NOX::Abstract::Group::ReturnType LOCA::BorderedSolver::Bordering::initForSolve ( )
virtual

Intialize solver for a solve.

This should be called after setMatrixBlocks(), but before applyInverse().

Implements LOCA::BorderedSolver::AbstractStrategy.

References NOX::Abstract::Group::Ok.

NOX::Abstract::Group::ReturnType LOCA::BorderedSolver::Bordering::initForTransposeSolve ( )
virtual

Intialize solver for a transpose solve.

This should be called after setMatrixBlocks(), but before applyInverseTranspose().

Implements LOCA::BorderedSolver::AbstractStrategy.

References NOX::Abstract::Group::Ok.

void LOCA::BorderedSolver::Bordering::setMatrixBlocks ( const Teuchos::RCP< const LOCA::BorderedSolver::AbstractOperator > &  op,
const Teuchos::RCP< const NOX::Abstract::MultiVector > &  blockA,
const Teuchos::RCP< const LOCA::MultiContinuation::ConstraintInterface > &  blockB,
const Teuchos::RCP< const NOX::Abstract::MultiVector::DenseMatrix > &  blockC 
)
virtual

Set blocks.

The blockA or blockC pointer may be null if either is zero. Whether block B is zero will be determined by querying blockB via ConstraintInterface::isConstraintDerivativesXZero.

Implements LOCA::BorderedSolver::AbstractStrategy.


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