| 
    NOX
    Development
    
   | 
 
Implementation of LOCA::MultiContinuation::ConstraintInterfaceMVDX for arclength continuation. More...
#include <LOCA_MultiContinuation_ArcLengthConstraint.H>
Public Member Functions | |
| ArcLengthConstraint (const Teuchos::RCP< LOCA::GlobalData > &global_data, const Teuchos::RCP< LOCA::MultiContinuation::ArcLengthGroup > &grp) | |
| Constructor.  | |
| ArcLengthConstraint (const ArcLengthConstraint &source, NOX::CopyType type=NOX::DeepCopy) | |
| Copy constructor.  | |
| ~ArcLengthConstraint () | |
| Destructor.  | |
| virtual void | setArcLengthGroup (const Teuchos::RCP< LOCA::MultiContinuation::ArcLengthGroup > &grp) | 
| Set pointer to arclength group.  | |
Implementation of LOCA::MultiContinuation::ConstraintInterface  | |
virtual methods  | |
| virtual void | copy (const 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.  | |
  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... | |
| virtual void | postProcessContinuationStep (LOCA::Abstract::Iterator::StepStatus) | 
| Perform any postprocessing after a continuation step finishes.  More... | |
Protected Attributes | |
| Teuchos::RCP< LOCA::GlobalData > | globalData | 
| LOCA global data object.  | |
| 
Teuchos::RCP < LOCA::MultiContinuation::ArcLengthGroup >  | arcLengthGroup | 
| Pointer to arc-length group.  | |
| NOX::Abstract::MultiVector::DenseMatrix | constraints | 
| Constraint values.  | |
| bool | isValidConstraints | 
| Flag indicating whether constraints are valid.  | |
| std::vector< int > | conParamIDs | 
| Continuation parameter IDs.  | |
Implementation of LOCA::MultiContinuation::ConstraintInterfaceMVDX for arclength continuation.
This class implements the arclength constraint equation for pseudo-arclength continuation:
 where 
, 
 are the solution and parameter components of the predictor direction 
 respectively. Since the derivative of 
 with respect to 
 is just 
, the predictor tangent, this class implements the MVDX version of the constraint interface. 
      
  | 
  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 LOCA::Extended::MultiVector::getScalar(), and NOX::Abstract::Group::Ok.
 1.8.5