|
NOX Development
|
Moore-Spence turning point solver strategy based on "Phipps" bordering which is the 5-solve modified turning point bordering algorithm that uses bordered linear solves. More...
#include <LOCA_TurningPoint_MooreSpence_PhippsBordering.H>


Public Member Functions | |
| PhippsBordering (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &solverParams) | |
| Constructor. | |
| virtual | ~PhippsBordering () |
| Destructor. | |
| virtual void | setBlocks (const Teuchos::RCP< LOCA::TurningPoint::MooreSpence::AbstractGroup > &group, const Teuchos::RCP< LOCA::TurningPoint::MooreSpence::ExtendedGroup > &tpGroup, const Teuchos::RCP< const NOX::Abstract::Vector > &nullVector, const Teuchos::RCP< const NOX::Abstract::Vector > &JnVector, const Teuchos::RCP< const NOX::Abstract::MultiVector > &dfdp, const Teuchos::RCP< const NOX::Abstract::MultiVector > &dJndp) |
| Set blocks in extended linear system. | |
| virtual NOX::Abstract::Group::ReturnType | solve (Teuchos::ParameterList ¶ms, const LOCA::TurningPoint::MooreSpence::ExtendedMultiVector &input, LOCA::TurningPoint::MooreSpence::ExtendedMultiVector &result) const |
| Solves the extended system as defined above. | |
| virtual NOX::Abstract::Group::ReturnType | solveTranspose (Teuchos::ParameterList ¶ms, const LOCA::TurningPoint::MooreSpence::ExtendedMultiVector &input, LOCA::TurningPoint::MooreSpence::ExtendedMultiVector &result) const |
| Solves the transpose of the extended system as defined above. | |
Public Member Functions inherited from LOCA::TurningPoint::MooreSpence::SolverStrategy | |
| SolverStrategy () | |
| Constructor. | |
| virtual | ~SolverStrategy () |
| Destructor. | |
Protected Member Functions | |
| NOX::Abstract::Group::ReturnType | solveContiguous (Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input_x, const NOX::Abstract::MultiVector &input_null, const NOX::Abstract::MultiVector::DenseMatrix &input_param, NOX::Abstract::MultiVector &result_x, NOX::Abstract::MultiVector &result_null, NOX::Abstract::MultiVector::DenseMatrix &result_param) const |
| Solves equations with contiguous arguments. | |
| NOX::Abstract::Group::ReturnType | solveTransposeContiguous (Teuchos::ParameterList ¶ms, const NOX::Abstract::MultiVector &input_x, const NOX::Abstract::MultiVector &input_null, const NOX::Abstract::MultiVector::DenseMatrix &input_param, NOX::Abstract::MultiVector &result_x, NOX::Abstract::MultiVector &result_null, NOX::Abstract::MultiVector::DenseMatrix &result_param) const |
| Solves equations with contiguous arguments. | |
Protected Attributes | |
| Teuchos::RCP< LOCA::GlobalData > | globalData |
| Global data object. | |
| Teuchos::RCP< Teuchos::ParameterList > | solverParams |
| Solver parameters. | |
| Teuchos::RCP< LOCA::TurningPoint::MooreSpence::AbstractGroup > | group |
| Underlying group. | |
| Teuchos::RCP< LOCA::TurningPoint::MooreSpence::ExtendedGroup > | tpGroup |
| Turning point group. | |
| Teuchos::RCP< const NOX::Abstract::Vector > | nullVector |
| Null vector. | |
| Teuchos::RCP< const NOX::Abstract::Vector > | JnVector |
| Jacobian times null vector. | |
| Teuchos::RCP< const NOX::Abstract::MultiVector > | dfdp |
| df/dp | |
| Teuchos::RCP< const NOX::Abstract::MultiVector > | dJndp |
| d(Jn)/dp | |
| Teuchos::RCP< LOCA::BorderedSolver::AbstractStrategy > | borderedSolver |
| Bordered solver for solving (n+1)x(n+1) sets of equations. | |
| Teuchos::RCP< LOCA::BorderedSolver::AbstractStrategy > | transposeBorderedSolver |
| Bordered solver for solving transpose (n+1)x(n+1) sets of equations. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | nullMultiVector |
| Null vector as a multivector. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | JnMultiVector |
| Jacobian times null vector as a multivector. | |
| Teuchos::RCP< const NOX::Abstract::Vector > | uVector |
| u vector | |
| Teuchos::RCP< const NOX::Abstract::MultiVector > | uMultiVector |
| u multi vector | |
| double | s |
| Norm of JnVector. | |
| double | st |
| Norm of J^T*u. | |
Moore-Spence turning point solver strategy based on "Phipps" bordering which is the 5-solve modified turning point bordering algorithm that uses bordered linear solves.
This class solves the Moore-Spence turning point Newton equations:
![\[
\begin{bmatrix}
J & 0 & f_p \\
(Jv)_x & J & (Jv)_p \\
0 & \phi^T & 0
\end{bmatrix}
\begin{bmatrix}
X \\
Y \\
z
\end{bmatrix} =
\begin{bmatrix}
F \\
G \\
h
\end{bmatrix}
\]](form_654.png)
via the following modified block elimination scheme:
![\[
\begin{split}
\begin{bmatrix}
J & u \\
v^T & 0
\end{bmatrix}
\begin{bmatrix}
A & B \\
a & b
\end{bmatrix} &=
\begin{bmatrix}
F & f_p \\
0 & 0
\end{bmatrix} \\
\begin{bmatrix}
J & u \\
v^T & 0
\end{bmatrix}
\begin{bmatrix}
C & D & E\\
c & d & e
\end{bmatrix} &=
\begin{bmatrix}
(Jv)_x A - G & (Jv)_x B - (Jv)_p & (Jv)_x v \\
0 & 0 & 0
\end{bmatrix} \\
\begin{bmatrix}
\sigma & 0 & b \\
e & \sigma & -d \\
-\phi^T E & \phi^T v & \phi^T D
\end{bmatrix}
\begin{bmatrix}
\alpha \\
\beta \\
z
\end{bmatrix} &=
\begin{bmatrix}
a \\
c \\
h + \phi^T C
\end{bmatrix} \\
X &= A - B z + v \alpha \\
Y &= -C + d z - E \alpha + v \beta
\end{split}
\]](form_655.png)
where 

| LOCA::TurningPoint::MooreSpence::PhippsBordering::PhippsBordering | ( | 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. Instantiates a bordered solver for solving the bordeded systems described above. See LOCA::BorderedSolver::Factory for a description of available solvers. |
References borderedSolver, globalData, and solverParams.
|
virtual |
Set blocks in extended linear system.
| group | [in] Underlying group representing J |
| tpGroup | [in] Turning point group representing the turning point equations. |
| nullVector | [in] Vector representing v |
| JnVector | [in] Vector representing Jv |
| dfdp | [in] Vector representing df/dp |
| dJndp | [in] Vector representing d(Jv)/dp |
Implements LOCA::TurningPoint::MooreSpence::SolverStrategy.
References NOX::DeepCopy, and NOX::Abstract::Vector::TwoNorm.
|
virtual |
Solves the extended system as defined above.
The params argument is the linear solver parameters.
Implements LOCA::TurningPoint::MooreSpence::SolverStrategy.
References LOCA::TurningPoint::MooreSpence::ExtendedMultiVector::getNullMultiVec(), LOCA::Extended::MultiVector::getScalars(), LOCA::TurningPoint::MooreSpence::ExtendedMultiVector::getXMultiVec(), and LOCA::Extended::MultiVector::numVectors().
|
virtual |
Solves the transpose of the extended system as defined above.
The params argument is the linear solver parameters.
Reimplemented from LOCA::TurningPoint::MooreSpence::SolverStrategy.
References LOCA::TurningPoint::MooreSpence::ExtendedMultiVector::getNullMultiVec(), LOCA::Extended::MultiVector::getScalars(), LOCA::TurningPoint::MooreSpence::ExtendedMultiVector::getXMultiVec(), and LOCA::Extended::MultiVector::numVectors().