|
NOX Development
|
Bordered system solver strategy based on bordering. More...
#include <LOCA_BorderedSolver_Bordering.H>


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. | |
| 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. | |
| virtual NOX::Abstract::Group::ReturnType | initForSolve () |
| Intialize solver for a solve. | |
| virtual NOX::Abstract::Group::ReturnType | initForTransposeSolve () |
| Intialize solver for a transpose solve. | |
| 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. | |
| 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. | |
| virtual NOX::Abstract::Group::ReturnType | applyInverse (Teuchos::ParameterList ¶ms, 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. | |
| virtual NOX::Abstract::Group::ReturnType | applyInverseTranspose (Teuchos::ParameterList ¶ms, 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. | |
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. | |
Protected Attributes | |
| Teuchos::RCP< LOCA::GlobalData > | globalData |
| 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 | |
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}
\]](form_314.png)
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}
\]](form_427.png)
It takes advantage of any of the matrix blocks being zero and concatenates 



To solve the transpose of the system, a similar bordering algorithm is implemented. Note however that for the transpose, the constraint object representing 

| 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.
| global_data | [in] Global data object |
| topParams | [in] Parsed top-level parameter list |
| solverParams | [in] Bordered solver parameters. Currently none are referenced. |
|
virtual |
Computed extended matrix-multivector product.
Computes
![\[
\begin{bmatrix}
U \\
V
\end{bmatrix} =
\begin{bmatrix}
J & A \\
B^T & C
\end{bmatrix}
\begin{bmatrix}
X \\
Y
\end{bmatrix} =
\begin{bmatrix}
J*X + A*Y \\
B^T*X + C*Y
\end{bmatrix}.
\]](form_316.png)
Implements LOCA::BorderedSolver::AbstractStrategy.
References NOX::Abstract::Group::Failed, NOX::Abstract::MultiVector::multiply(), and NOX::Abstract::MultiVector::update().
|
virtual |
Solves the extended system as defined above using bordering.
The params argument is the linear solver parameters. If isZeroF or isZeroG is true, than the corresponding F or G pointers may be NULL.
Implements LOCA::BorderedSolver::AbstractStrategy.
References NOX::Abstract::MultiVector::clone(), NOX::Abstract::MultiVector::numVectors(), LOCA::BorderedSolver::LowerTriangularBlockElimination::solve(), and LOCA::BorderedSolver::UpperTriangularBlockElimination::solve().
|
virtual |
Solves the transpose of the extended system as defined above using bordering.
The params argument is the linear solver parameters. If isZeroF or isZeroG is true, than the corresponding F or G pointers may be NULL. Note that for the transpose solve B must be of type LOCA::MultiContinuation::ConstraintInterfaceMVDX.
Implements LOCA::BorderedSolver::AbstractStrategy.
References NOX::Abstract::MultiVector::clone(), NOX::Abstract::MultiVector::numVectors(), LOCA::BorderedSolver::LowerTriangularBlockElimination::solveTranspose(), and LOCA::BorderedSolver::UpperTriangularBlockElimination::solveTranspose().
|
virtual |
Computed extended matrix transpose-multivector product.
Computes
![\[
\begin{bmatrix}
U \\
V
\end{bmatrix} =
\begin{bmatrix}
J^T & B \\
A^T & C
\end{bmatrix}
\begin{bmatrix}
X \\
Y
\end{bmatrix} =
\begin{bmatrix}
J^T*X + B*Y \\
A^T*X + C^T*Y
\end{bmatrix}.
\]](form_317.png)
Implements LOCA::BorderedSolver::AbstractStrategy.
References NOX::Abstract::Group::Failed, and NOX::Abstract::MultiVector::multiply().
|
virtual |
Intialize solver for a solve.
This should be called after setMatrixBlocks(), but before applyInverse().
Implements LOCA::BorderedSolver::AbstractStrategy.
References NOX::Abstract::Group::Ok.
|
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.
|
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.