|
NOX Development
|
Implementation of LOCA::MultiContinuation::ConstraintInterfaceMVDX for computing turning points for the minimally augmented turning point formulation. More...
#include <LOCA_TurningPoint_MinimallyAugmented_ModifiedConstraint.H>


Public Member Functions | |
| ModifiedConstraint (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &tpParams, const Teuchos::RCP< LOCA::TurningPoint::MinimallyAugmented::AbstractGroup > &g, int bif_param) | |
| Constructor. | |
| ModifiedConstraint (const ModifiedConstraint &source, NOX::CopyType type=NOX::DeepCopy) | |
| Copy constructor. | |
| virtual | ~ModifiedConstraint () |
| Destructor. | |
Public Member Functions inherited from LOCA::TurningPoint::MinimallyAugmented::Constraint | |
| Constraint (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &tpParams, const Teuchos::RCP< LOCA::TurningPoint::MinimallyAugmented::AbstractGroup > &g, int bif_param) | |
| Constructor. | |
| Constraint (const Constraint &source, NOX::CopyType type=NOX::DeepCopy) | |
| Copy constructor. | |
| virtual | ~Constraint () |
| Destructor. | |
| virtual void | setGroup (const Teuchos::RCP< LOCA::TurningPoint::MinimallyAugmented::AbstractGroup > &g) |
| Set the group pointer. | |
| virtual Teuchos::RCP< const NOX::Abstract::Vector > | getLeftNullVec () const |
| Returns left null vector w. | |
| virtual Teuchos::RCP< const NOX::Abstract::Vector > | getRightNullVec () const |
| Returns right null vector v. | |
| virtual Teuchos::RCP< const NOX::Abstract::Vector > | getAVec () const |
| Returns a vector. | |
| virtual Teuchos::RCP< const NOX::Abstract::Vector > | getBVec () const |
| Returns b vector. | |
| virtual double | getSigma () const |
| Returns sigma. | |
| virtual int | numConstraints () const |
| Return number of constraints. | |
| virtual void | setX (const NOX::Abstract::Vector &y) |
| Set the solution vector to y. | |
| virtual void | setParam (int paramID, double val) |
| Sets parameter indexed by paramID. | |
| virtual void | setParams (const std::vector< int > ¶mIDs, const NOX::Abstract::MultiVector::DenseMatrix &vals) |
| Sets parameters indexed by paramIDs. | |
| virtual NOX::Abstract::Group::ReturnType | computeDX () |
| Compute derivative of constraints w.r.t. solution vector x. | |
| virtual NOX::Abstract::Group::ReturnType | computeDP (const std::vector< int > ¶mIDs, NOX::Abstract::MultiVector::DenseMatrix &dgdp, bool isValidG) |
| Compute derivative of constraints w.r.t. supplied parameters. | |
| virtual bool | isConstraints () const |
Return true if constraint residuals are valid. | |
| virtual bool | isDX () const |
Return true if derivatives of constraints w.r.t. x are valid. | |
| virtual const NOX::Abstract::MultiVector::DenseMatrix & | getConstraints () const |
| Return constraint residuals. | |
| virtual const NOX::Abstract::MultiVector * | getDX () const |
| Return solution component of constraint derivatives. | |
| virtual bool | isDXZero () const |
Return true if solution component of constraint derivatives is zero. | |
Public Member Functions inherited from LOCA::MultiContinuation::ConstraintInterfaceMVDX | |
| ConstraintInterfaceMVDX () | |
| Constructor. | |
| virtual | ~ConstraintInterfaceMVDX () |
| Destructor. | |
| virtual NOX::Abstract::Group::ReturnType | multiplyDX (double alpha, const NOX::Abstract::MultiVector &input_x, NOX::Abstract::MultiVector::DenseMatrix &result_p) const |
| Compute result_p = alpha * dg/dx * input_x. | |
| virtual NOX::Abstract::Group::ReturnType | addDX (Teuchos::ETransp transb, double alpha, const NOX::Abstract::MultiVector::DenseMatrix &b, double beta, NOX::Abstract::MultiVector &result_x) const |
| Compute result_x = alpha * dg/dx^T * op(b) + beta * result_x. | |
Public Member Functions inherited from LOCA::MultiContinuation::ConstraintInterface | |
| ConstraintInterface () | |
| Constructor. | |
| virtual | ~ConstraintInterface () |
| Destructor. | |
Implementation of LOCA::MultiContinuation::ConstraintInterface | |
virtual methods | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | w_vector_update |
| Stores update to left null vector. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | v_vector_update |
| Stores update to right null vector. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | w_residual |
| Stores left null vector residual. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | v_residual |
| Stores right null vector residual. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | deltaX |
| Stores solution update. | |
| NOX::Abstract::MultiVector::DenseMatrix | sigma1 |
| Stores sigma_1. | |
| NOX::Abstract::MultiVector::DenseMatrix | sigma2 |
| Stores sigma_1. | |
| double | deltaP |
| Stores parameter update. | |
| bool | isFirstSolve |
| bool | includeNewtonTerms |
| Flag indicating whether to include the newton update terms. | |
| virtual void | copy (const LOCA::MultiContinuation::ConstraintInterface &source) |
| Copy. | |
| virtual Teuchos::RCP< LOCA::MultiContinuation::ConstraintInterface > | clone (NOX::CopyType type=NOX::DeepCopy) const |
| Cloning function. | |
| virtual NOX::Abstract::Group::ReturnType | computeConstraints () |
| Compute continuation constraint equations. | |
| virtual void | preProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus stepStatus) |
| Perform any preprocessing before a continuation step starts. | |
| virtual void | postProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus stepStatus) |
| Perform any postprocessing after a continuation step finishes. | |
| void | setNewtonUpdates (const NOX::Abstract::Vector &dx, double dp, double step) |
| Set the newton update for x and p. | |
Additional Inherited Members | |
Protected Types inherited from LOCA::TurningPoint::MinimallyAugmented::Constraint | |
| enum | NullVectorScaling { NVS_None , NVS_OrderOne , NVS_OrderN } |
| Enumerated type determining type of scaling. More... | |
Protected Member Functions inherited from LOCA::TurningPoint::MinimallyAugmented::Constraint | |
| virtual void | scaleNullVectors (NOX::Abstract::Vector &a, NOX::Abstract::Vector &b) |
| Scale a & b vectors. | |
| virtual void | getInitialVectors (NOX::Abstract::Vector &a, NOX::Abstract::Vector &b) |
| Get initial a & b vectors. | |
Protected Attributes inherited from LOCA::TurningPoint::MinimallyAugmented::Constraint | |
| Teuchos::RCP< LOCA::GlobalData > | globalData |
| Pointer LOCA global data object. | |
| Teuchos::RCP< LOCA::Parameter::SublistParser > | parsedParams |
| Parsed top-level parameters. | |
| Teuchos::RCP< Teuchos::ParameterList > | turningPointParams |
| Bifurcation parameter list. | |
| Teuchos::RCP< LOCA::TurningPoint::MinimallyAugmented::AbstractGroup > | grpPtr |
| Pointer to base group that defines | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | a_vector |
| Vector for | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | b_vector |
| Vector for | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | w_vector |
| Stores left null vector. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | v_vector |
| Stores right null vector. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | Jv_vector |
| Stores J*v. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | Jtw_vector |
| Stores J^T*w. | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | sigma_x |
| Stores sigma_x. | |
| NOX::Abstract::MultiVector::DenseMatrix | constraints |
| Constraint values. | |
| Teuchos::RCP< LOCA::BorderedSolver::AbstractStrategy > | borderedSolver |
| Stores bordered solver strategy. | |
| double | dn |
| Stores vector length as a double. | |
| double | sigma_scale |
| Stores scale factor on sigma. | |
| bool | isSymmetric |
| Flag indicating whether Jacobian is symmetric. | |
| bool | isValidConstraints |
| Flag indicating whether constraints are valid. | |
| bool | isValidDX |
| Flag indicating whether sigma_x is valid. | |
| std::vector< int > | bifParamID |
| Stores the bifurcation parameter index. | |
| bool | updateVectorsEveryContinuationStep |
| Flag indicating whether to update | |
| bool | updateVectorsEveryIteration |
| Flag indicating whether to update | |
| NullVectorScaling | nullVecScaling |
| Null vector scaling method. | |
| bool | multiplyMass |
| Multiply null vectors by mass matrix. | |
| Teuchos::RCP< LOCA::TimeDependent::AbstractGroup > | tdGrp |
| Time dependent group interface for multiplying by mass matrix. | |
| Teuchos::RCP< NOX::Abstract::Vector > | tmp_mass |
| Temporary vector for mulplying null vectors by mass matrix. | |
Implementation of LOCA::MultiContinuation::ConstraintInterfaceMVDX for computing turning points for the minimally augmented turning point formulation.
This class is a modification of LOCA::TurningPoint::MinimallyAugmented::Constraint where updates are computed to the left and right null vectors 

![\[
\begin{bmatrix}
J & a \\
b^T & 0
\end{bmatrix}
\begin{bmatrix}
\Delta v \\
\Delta \sigma_1
\end{bmatrix} = -
\begin{bmatrix}
J v + \sigma_1 a + (Jv)_x\Delta x + (Jv)_p \Delta p\\
b^T a - n
\end{bmatrix},
\]](form_642.png)
![\[
\begin{bmatrix}
J^T & b \\
a^T & 0
\end{bmatrix}
\begin{bmatrix}
\Delta w \\
\Delta \sigma_2
\end{bmatrix} = -
\begin{bmatrix}
J^T w + \sigma_2 b + (J^T w)_x \Delta x + (J^T w)_p \Delta p \\
a^T w - n
\end{bmatrix},
\]](form_643.png)
The class is intialized via the tpParams parameter list argument to the constructor. This class recognizes all paramters for LOCA::TurningPoint::MinimallyAugmented::Constraint plus the following:


|
virtual |
Cloning function.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
|
virtual |
Compute continuation constraint equations.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
References NOX::Abstract::Group::Ok, NOX::Utils::OuterIteration, and NOX::ShapeCopy.
|
virtual |
Copy.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
References LOCA::TurningPoint::MinimallyAugmented::Constraint::copy(), deltaP, deltaX, includeNewtonTerms, isFirstSolve, sigma1, sigma2, v_residual, v_vector_update, w_residual, and w_vector_update.
|
virtual |
Perform any postprocessing after a continuation step finishes.
The stepStatus argument indicates whether the step was successful. Here we set up the constraint class to solve for 

Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
References LOCA::TurningPoint::MinimallyAugmented::Constraint::postProcessContinuationStep(), and LOCA::Abstract::Iterator::Successful.
|
virtual |
Perform any preprocessing before a continuation step starts.
The stepStatus argument indicates whether the previous step was successful. Here we set up the constraint class to solve for 

Reimplemented from LOCA::MultiContinuation::ConstraintInterface.
References LOCA::MultiContinuation::ConstraintInterface::preProcessContinuationStep().
|
protected |
Flag that indicates whether we're in the first solve per continuation step.
Referenced by copy().