NOX
Development
|
Implementation of LOCA::MultiContinuation::ConstraintInterfaceMVDX for computing Hopf bifurcations for the minimally augmented Hopf formulation. More...
#include <LOCA_Hopf_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 > &hpfParams, const Teuchos::RCP< LOCA::Hopf::MinimallyAugmented::AbstractGroup > &g, bool is_symmetric, const NOX::Abstract::Vector &a_real, const NOX::Abstract::Vector &a_imag, const NOX::Abstract::Vector *b_real, const NOX::Abstract::Vector *b_imag, int bif_param, double freq) | |
Constructor. | |
Constraint (const Constraint &source, NOX::CopyType type=NOX::DeepCopy) | |
Copy constructor. | |
virtual | ~Constraint () |
Destructor. | |
virtual void | setGroup (const Teuchos::RCP< LOCA::Hopf::MinimallyAugmented::AbstractGroup > &g) |
Set the group pointer. More... | |
virtual void | setFrequency (double freq) |
Set Hopf frequency. | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getLeftNullVecReal () const |
Returns real component of left null vector w. | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getLeftNullVecImag () const |
Returns imaginary component of left null vector w. | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getRightNullVecReal () const |
Returns real component of right null vector v. | |
virtual Teuchos::RCP< const NOX::Abstract::Vector > | getRightNullVecImag () const |
Returns imaginary component of right null vector v. | |
virtual double | getSigmaReal () const |
Returns real component of sigma. | |
virtual double | getSigmaImag () const |
Returns imaginary component of sigma. | |
virtual NOX::Abstract::Group::ReturnType | computeDOmega (NOX::Abstract::MultiVector::DenseMatrix &domega) |
Compute derivative of sigma w.r.t. frequency omega. | |
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 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. | |
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 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 > | hopfParams |
Bifurcation parameter list. | |
Teuchos::RCP < LOCA::Hopf::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 > | Cv_vector |
Stores C*v. | |
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. | |
double | omega |
Stores Hopf frequency. | |
bool | updateVectorsEveryContinuationStep |
Flag indicating whether to update and every continuation step. | |
bool | updateVectorsEveryIteration |
Flag indicating whether to update and every nonlinear iteration. | |
Implementation of LOCA::MultiContinuation::ConstraintInterfaceMVDX for computing Hopf bifurcations for the minimally augmented Hopf formulation.
This class implements the turning point constraint equation for the minimally augmented Hopf formulation where is defined via
for any vectors and in . Using these relationships, it is easy to show
The class is intialized via the hpfParams
parameter list argument to the constructor. The parameters this class recognizes are:
|
virtual |
Compute derivative of constraints w.r.t. supplied parameters.
The first column of dgdp
should be filled with the constraint residuals if isValidG
is false
. If isValidG
is true
, then the dgdp
contains on input.
Implements LOCA::MultiContinuation::ConstraintInterface.
References NOX::Abstract::Group::Ok, Teuchos::SerialDenseMatrix< OrdinalType, ScalarType >::scale(), and Teuchos::View.
|
virtual |
Perform any postprocessing after a continuation step finishes.
The stepStatus
argument indicates whether the step was successful. Here we update the and vectors to and respectively if requested.
Reimplemented from LOCA::MultiContinuation::ConstraintInterface.
References NOX::Utils::StepperDetails, and LOCA::Abstract::Iterator::Successful.
|
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.