| NOX
    Development
    | 
Implementation of LOCA::MultiContinuation::ConstraintInterfaceMVDX for computing turning points for the minimally augmented turning point formulation. More...
#include <LOCA_TurningPoint_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 > &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.  More... | |
| 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. | |
| Implementation of LOCA::MultiContinuation::ConstraintInterface | |
| virtual methods | |
| 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 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 | 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.  More... | |
| virtual bool | isConstraints () const | 
| Return trueif constraint residuals are valid. | |
| virtual bool | isDX () const | 
| Return trueif 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 trueif solution component of constraint derivatives is zero. | |
| virtual void | postProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus stepStatus) | 
| Perform any postprocessing after a continuation step finishes.  More... | |
|  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.  More... | |
| 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.  More... | |
|  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.  More... | |
| Protected Types | |
| enum | NullVectorScaling { NVS_None, NVS_OrderOne, NVS_OrderN } | 
| Enumerated type determining type of scaling.  More... | |
| Protected Member Functions | |
| 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 | |
| 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  and  every continuation step. | |
| bool | updateVectorsEveryIteration | 
| Flag indicating whether to update  and  every nonlinear iteration. | |
| 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 implements the turning point constraint equation  for the minimally augmented turning point formulation where
 for the minimally augmented turning point formulation where  is defined via
 is defined via 
![\[ \begin{bmatrix} J & a \\ b^T & 0 \end{bmatrix} \begin{bmatrix} v \\ \sigma_1 \end{bmatrix} = \begin{bmatrix} 0 \\ n \end{bmatrix}, \]](form_503.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_504.png) 
![\[ \sigma = -w^T J v/n \]](form_505.png) 
 for any vectors  and
 and  in
 in  . Using these relationships, it is easy to show
. Using these relationships, it is easy to show 
![\[ \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_507.png) 
The class is intialized via the tpParams parameter list argument to the constructor. The parameters this class recognizes are: 
 is symmetric, in which case we force
 is symmetric, in which case we force  and therefore the second tranpose solve for
 and therefore the second tranpose solve for  is unnecessary
 is unnecessary  and
 and  vectors. Valid choices are:
 vectors. Valid choices are:  vector
 vector  vector
 vector  and
 and  where
 where  is the bifurcation parameter.
 is the bifurcation parameter.  and
 and  are set to 1.0
 are set to 1.0  and
 and  . This determines the norm of these vectors and the scaling of
. This determines the norm of these vectors and the scaling of  . Valid choices are:
. Valid choices are:  and
 and  vectors via
 vectors via  and
 and  every continuation step
 every continuation step  and
 and  vectors via
 vectors via  and
 and  every nonlinear iteration
 every nonlinear iteration  and
 and  vectors by the mass matrix
 vectors by the mass matrix  at the strart of a turning point calculation, and each time
 at the strart of a turning point calculation, and each time  and
 and  are updated. This can improve the scaling of these vectors, and may orthogonalize them against structural null spaces (i.e., pressure null space for incompressible Navier-Stokes).
 are updated. This can improve the scaling of these vectors, and may orthogonalize them against structural null spaces (i.e., pressure null space for incompressible Navier-Stokes). | 
 | virtual | 
Compute derivative of constraints w.r.t. supplied parameters.
The first column of dgdp should be filled with the constraint residuals  if
 if isValidG is false. If isValidG is true, then the dgdp contains  on input.
 on input. 
Implements LOCA::MultiContinuation::ConstraintInterface.
Reimplemented in LOCA::Pitchfork::MinimallyAugmented::Constraint.
References NOX::Abstract::Group::Ok, and Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::scale().
Referenced by LOCA::Pitchfork::MinimallyAugmented::Constraint::computeDP().
| 
 | virtual | 
Perform any postprocessing after a continuation step finishes.
The stepStatus argument indicates whether the step was successful. Here we update the  and
 and  vectors to
 vectors to  and
 and  respectively if requested.
 respectively if requested. 
Reimplemented from LOCA::MultiContinuation::ConstraintInterface.
Reimplemented in LOCA::TurningPoint::MinimallyAugmented::ModifiedConstraint.
References NOX::Utils::StepperDetails, and LOCA::Abstract::Iterator::Successful.
Referenced by LOCA::TurningPoint::MinimallyAugmented::ModifiedConstraint::postProcessContinuationStep().
| 
 | 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.
Reimplemented in LOCA::Pitchfork::MinimallyAugmented::Constraint.
Referenced by LOCA::TurningPoint::MinimallyAugmented::ExtendedGroup::copy(), LOCA::TurningPoint::MinimallyAugmented::ExtendedGroup::ExtendedGroup(), and LOCA::Pitchfork::MinimallyAugmented::Constraint::setGroup().
 1.8.5
 1.8.5