|
NOX Development
|
Implementation of LOCA::MultiContinuation::ConstraintInterfaceMVDX for computing pitchforks for the minimally augmented pitchfork formulation. More...
#include <LOCA_Pitchfork_MinimallyAugmented_Constraint.H>


Public Member Functions | |
| Constraint (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::Parameter::SublistParser > &topParams, const Teuchos::RCP< Teuchos::ParameterList > &pfParams, const Teuchos::RCP< LOCA::Pitchfork::MinimallyAugmented::AbstractGroup > &g, const Teuchos::RCP< const NOX::Abstract::Vector > &psi, 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. | |
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 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 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 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 bool | isDXZero () const |
Return true if solution component of constraint derivatives is zero. | |
| virtual void | postProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus stepStatus) |
| Perform any postprocessing after a continuation step finishes. | |
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. | |
| virtual void | preProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus) |
| Perform any preprocessing before a continuation step starts. | |
Implementation of LOCA::MultiContinuation::ConstraintInterface | |
virtual methods | |
| Teuchos::RCP< LOCA::Pitchfork::MinimallyAugmented::AbstractGroup > | pf_grp |
| Pointer to pitchfork group. | |
| Teuchos::RCP< const NOX::Abstract::Vector > | psi_vector |
| Vector for | |
| Teuchos::RCP< NOX::Abstract::MultiVector > | dgdx |
| Stores dg/dx. | |
| NOX::Abstract::MultiVector::DenseMatrix | pf_constraints |
| Constraint values. | |
| 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 int | numConstraints () const |
| Return number of constraints. | |
| virtual NOX::Abstract::Group::ReturnType | computeConstraints () |
| Compute continuation constraint equations. | |
| 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 const NOX::Abstract::MultiVector::DenseMatrix & | getConstraints () const |
| Return constraint residuals. | |
| virtual const NOX::Abstract::MultiVector * | getDX () const |
| Return solution component of constraint derivatives. | |
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 pitchforks for the minimally augmented pitchfork formulation.
This class implements the pitchfork constraint equations 


![\[
\begin{bmatrix}
J & a \\
b^T & 0
\end{bmatrix}
\begin{bmatrix}
v \\
\sigma_1
\end{bmatrix} =
\begin{bmatrix}
0 \\
n
\end{bmatrix},
\]](form_580.png)
![\[
\begin{bmatrix}
J^T & b \\
a^T & 0
\end{bmatrix}
\begin{bmatrix}
w \\
\sigma_2
\end{bmatrix} =
\begin{bmatrix}
0 \\
n
\end{bmatrix},
\]](form_581.png)
![\[
\sigma = -w^T J v/n
\]](form_582.png)
for any vectors 


![\[
\begin{split}
\sigma_x &= -(w^T J v)_x/n = -w^T J_x v/n \\
\sigma_p &= -(w^T J v)_p/n = -w^T J_p v/n
\end{split}
\]](form_584.png)
The class is derived from LOCA::TurningPoint::MinimallyAugmented::Constraint. See this class for a description of available parameters.
|
virtual |
Destructor.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
|
virtual |
Cloning function.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
|
virtual |
Compute continuation constraint equations.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
References LOCA::TurningPoint::MinimallyAugmented::Constraint::computeConstraints(), and NOX::Abstract::Group::Ok.
|
virtual |
Compute derivative of constraints w.r.t. supplied parameters.
The first column of dgdp should be filled with the constraint residuals 
isValidG is false. If isValidG is true, then the dgdp contains 
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
References LOCA::TurningPoint::MinimallyAugmented::Constraint::computeDP().
|
virtual |
Compute derivative of constraints w.r.t. solution vector x.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
References LOCA::TurningPoint::MinimallyAugmented::Constraint::computeDX(), and NOX::Abstract::Group::Ok.
|
virtual |
Copy.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
References LOCA::TurningPoint::MinimallyAugmented::Constraint::copy(), dgdx, pf_constraints, and psi_vector.
|
virtual |
Return constraint residuals.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
|
virtual |
Return solution component of constraint derivatives.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
|
virtual |
Return number of constraints.
Reimplemented from LOCA::TurningPoint::MinimallyAugmented::Constraint.
|
virtual |
Set the group pointer.
This method should be called when ever the constrained group is copied, since we don't explicitly copy the underlying group here.
References LOCA::TurningPoint::MinimallyAugmented::Constraint::setGroup().